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