scalarMatricesTemplates.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 -------------------------------------------------------------------------------
10 License
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 "ListOps.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class Type>
34 void Foam::solve
35 (
36  scalarSquareMatrix& tmpMatrix,
37  List<Type>& sourceSol
38 )
39 {
40  label m = tmpMatrix.m();
41 
42  // Elimination
43  for (label i = 0; i < m; ++i)
44  {
45  label iMax = i;
46  scalar largestCoeff = mag(tmpMatrix(iMax, i));
47 
48  // Swap elements around to find a good pivot
49  for (label j = i + 1; j < m; ++j)
50  {
51  if (mag(tmpMatrix(j, i)) > largestCoeff)
52  {
53  iMax = j;
54  largestCoeff = mag(tmpMatrix(iMax, i));
55  }
56  }
57 
58  if (i != iMax)
59  {
60  for (label k = i; k < m; ++k)
61  {
62  Foam::Swap(tmpMatrix(i, k), tmpMatrix(iMax, k));
63  }
64  Foam::Swap(sourceSol[i], sourceSol[iMax]);
65  }
66 
67  // Check that the system of equations isn't singular
68  if (mag(tmpMatrix(i, i)) < 1e-20)
69  {
71  << "Singular Matrix"
72  << exit(FatalError);
73  }
74 
75  // Reduce to upper triangular form
76  for (label j = i + 1; j < m; ++j)
77  {
78  sourceSol[j] -= sourceSol[i]*(tmpMatrix(j, i)/tmpMatrix(i, i));
79 
80  for (label k = m - 1; k >= i; --k)
81  {
82  tmpMatrix(j, k) -=
83  tmpMatrix(i, k)*tmpMatrix(j, i)/tmpMatrix(i, i);
84  }
85  }
86  }
87 
88  // Back-substitution
89  for (label j = m - 1; j >= 0; --j)
90  {
91  Type ntempvec = Zero;
92 
93  for (label k = j + 1; k < m; ++k)
94  {
95  ntempvec += tmpMatrix(j, k)*sourceSol[k];
96  }
97 
98  sourceSol[j] = (sourceSol[j] - ntempvec)/tmpMatrix(j, j);
99  }
100 }
101 
102 
103 template<class Type>
104 void Foam::solve
105 (
106  List<Type>& psi,
107  const scalarSquareMatrix& matrix,
108  const List<Type>& source
109 )
110 {
111  scalarSquareMatrix tmpMatrix = matrix;
112  psi = source;
113  solve(tmpMatrix, psi);
114 }
115 
116 
117 template<class Type>
119 (
120  const scalarSquareMatrix& luMatrix,
121  const labelList& pivotIndices,
122  List<Type>& sourceSol
123 )
124 {
125  label m = luMatrix.m();
126 
127  label ii = 0;
128 
129  for (label i = 0; i < m; ++i)
130  {
131  label ip = pivotIndices[i];
132  Type sum = sourceSol[ip];
133  sourceSol[ip] = sourceSol[i];
134  const scalar* __restrict__ luMatrixi = luMatrix[i];
135 
136  if (ii != 0)
137  {
138  for (label j = ii - 1; j < i; ++j)
139  {
140  sum -= luMatrixi[j]*sourceSol[j];
141  }
142  }
143  else if (sum != Type(Zero))
144  {
145  ii = i + 1;
146  }
147 
148  sourceSol[i] = sum;
149  }
150 
151  for (label i = m - 1; i >= 0; --i)
152  {
153  Type sum = sourceSol[i];
154  const scalar* __restrict__ luMatrixi = luMatrix[i];
155 
156  for (label j = i + 1; j < m; ++j)
157  {
158  sum -= luMatrixi[j]*sourceSol[j];
159  }
160 
161  sourceSol[i] = sum/luMatrixi[i];
162  }
163 }
164 
165 
166 template<class Type>
168 (
169  const scalarSymmetricSquareMatrix& luMatrix,
170  List<Type>& sourceSol
171 )
172 {
173  label m = luMatrix.m();
174 
175  label ii = 0;
176 
177  for (label i = 0; i < m; ++i)
178  {
179  Type sum = sourceSol[i];
180  const scalar* __restrict__ luMatrixi = luMatrix[i];
181 
182  if (ii != 0)
183  {
184  for (label j = ii - 1; j < i; ++j)
185  {
186  sum -= luMatrixi[j]*sourceSol[j];
187  }
188  }
189  else if (sum != Type(Zero))
190  {
191  ii = i + 1;
192  }
193 
194  sourceSol[i] = sum/luMatrixi[i];
195  }
196 
197  for (label i = m - 1; i >= 0; --i)
198  {
199  Type sum = sourceSol[i];
200  const scalar* __restrict__ luMatrixi = luMatrix[i];
201 
202  for (label j = i + 1; j < m; ++j)
203  {
204  sum -= luMatrixi[j]*sourceSol[j];
205  }
206 
207  sourceSol[i] = sum/luMatrixi[i];
208  }
209 }
210 
211 
212 template<class Type>
213 void Foam::LUsolve
214 (
215  scalarSquareMatrix& matrix,
216  List<Type>& sourceSol
217 )
218 {
219  labelList pivotIndices(matrix.m());
220  LUDecompose(matrix, pivotIndices);
221  LUBacksubstitute(matrix, pivotIndices, sourceSol);
222 }
223 
224 
225 template<class Type>
226 void Foam::LUsolve
227 (
229  List<Type>& sourceSol
230 )
231 {
232  LUDecompose(matrix);
233  LUBacksubstitute(matrix, sourceSol);
234 }
235 
236 
237 template<class Form, class Type>
238 void Foam::multiply
239 (
240  Matrix<Form, Type>& ans, // value changed in return
241  const Matrix<Form, Type>& A,
242  const Matrix<Form, Type>& B
243 )
244 {
245  if (A.n() != B.m())
246  {
248  << "A and B must have identical inner dimensions but A.n = "
249  << A.n() << " and B.m = " << B.m()
250  << abort(FatalError);
251  }
252 
253  ans = Matrix<Form, Type>(A.m(), B.n(), Zero);
254 
255  for (label i = 0; i < A.m(); ++i)
256  {
257  for (label j = 0; j < B.n(); ++j)
258  {
259  for (label l = 0; l < B.m(); ++l)
260  {
261  ans(i, j) += A(i, l)*B(l, j);
262  }
263  }
264  }
265 }
266 
267 
268 // ************************************************************************* //
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:598
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
label k
Boltzmann constant.
Various functions to operate on Lists.
CEqn solve()
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
errorManip< error > abort(error &err)
Definition: errorManip.H:139
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
SymmetricSquareMatrix< scalar > scalarSymmetricSquareMatrix
label m() const noexcept
The number of rows.
Definition: Matrix.H:248
void Swap(DynamicList< T, SizeMinA > &a, DynamicList< T, SizeMinB > &b)
Exchange contents of lists - see DynamicList::swap().
Definition: DynamicList.H:692
const volScalarField & psi
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)
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
SquareMatrix< scalar > scalarSquareMatrix
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting returning the LU form of the matrix and the sol...
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution in the source.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127