DILUPreconditioner.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "DILUPreconditioner.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(DILUPreconditioner, 0);
33 
34  lduMatrix::preconditioner::
35  addasymMatrixConstructorToTable<DILUPreconditioner>
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const lduMatrix::solver& sol,
45  const dictionary&
46 )
47 :
49  rD_(sol.matrix().diag())
50 {
51  calcReciprocalD(rD_, sol.matrix());
52 }
53 
54 
55 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56 
58 (
59  scalarField& rD,
60  const lduMatrix& matrix
61 )
62 {
63  scalar* __restrict__ rDPtr = rD.begin();
64 
65  const label* const __restrict__ uPtr = matrix.lduAddr().upperAddr().begin();
66  const label* const __restrict__ lPtr = matrix.lduAddr().lowerAddr().begin();
67 
68  const scalar* const __restrict__ upperPtr = matrix.upper().begin();
69  const scalar* const __restrict__ lowerPtr = matrix.lower().begin();
70 
71  label nFaces = matrix.upper().size();
72  for (label face=0; face<nFaces; face++)
73  {
74  rDPtr[uPtr[face]] -= upperPtr[face]*lowerPtr[face]/rDPtr[lPtr[face]];
75  }
76 
77 
78  // Calculate the reciprocal of the preconditioned diagonal
79  label nCells = rD.size();
80 
81  for (label cell=0; cell<nCells; cell++)
82  {
83  rDPtr[cell] = 1.0/rDPtr[cell];
84  }
85 }
86 
87 
89 (
90  scalarField& wA,
91  const scalarField& rA,
92  const direction
93 ) const
94 {
95  scalar* __restrict__ wAPtr = wA.begin();
96  const scalar* __restrict__ rAPtr = rA.begin();
97  const scalar* __restrict__ rDPtr = rD_.begin();
98 
99  const label* const __restrict__ uPtr =
100  solver_.matrix().lduAddr().upperAddr().begin();
101  const label* const __restrict__ lPtr =
102  solver_.matrix().lduAddr().lowerAddr().begin();
103  const label* const __restrict__ losortPtr =
104  solver_.matrix().lduAddr().losortAddr().begin();
105 
106  const scalar* const __restrict__ upperPtr =
107  solver_.matrix().upper().begin();
108  const scalar* const __restrict__ lowerPtr =
109  solver_.matrix().lower().begin();
110 
111  label nCells = wA.size();
112  label nFaces = solver_.matrix().upper().size();
113  label nFacesM1 = nFaces - 1;
114 
115  for (label cell=0; cell<nCells; cell++)
116  {
117  wAPtr[cell] = rDPtr[cell]*rAPtr[cell];
118  }
119 
120 
121  label sface;
122 
123  for (label face=0; face<nFaces; face++)
124  {
125  sface = losortPtr[face];
126  wAPtr[uPtr[sface]] -=
127  rDPtr[uPtr[sface]]*lowerPtr[sface]*wAPtr[lPtr[sface]];
128  }
129 
130  for (label face=nFacesM1; face>=0; face--)
131  {
132  wAPtr[lPtr[face]] -=
133  rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[uPtr[face]];
134  }
135 }
136 
137 
139 (
140  scalarField& wT,
141  const scalarField& rT,
142  const direction
143 ) const
144 {
145  scalar* __restrict__ wTPtr = wT.begin();
146  const scalar* __restrict__ rTPtr = rT.begin();
147  const scalar* __restrict__ rDPtr = rD_.begin();
148 
149  const label* const __restrict__ uPtr =
150  solver_.matrix().lduAddr().upperAddr().begin();
151  const label* const __restrict__ lPtr =
152  solver_.matrix().lduAddr().lowerAddr().begin();
153  const label* const __restrict__ losortPtr =
154  solver_.matrix().lduAddr().losortAddr().begin();
155 
156  const scalar* const __restrict__ upperPtr =
157  solver_.matrix().upper().begin();
158  const scalar* const __restrict__ lowerPtr =
159  solver_.matrix().lower().begin();
160 
161  label nCells = wT.size();
162  label nFaces = solver_.matrix().upper().size();
163  label nFacesM1 = nFaces - 1;
164 
165  for (label cell=0; cell<nCells; cell++)
166  {
167  wTPtr[cell] = rDPtr[cell]*rTPtr[cell];
168  }
169 
170  for (label face=0; face<nFaces; face++)
171  {
172  wTPtr[uPtr[face]] -=
173  rDPtr[uPtr[face]]*upperPtr[face]*wTPtr[lPtr[face]];
174  }
175 
176 
177  label sface;
178 
179  for (label face=nFacesM1; face>=0; face--)
180  {
181  sface = losortPtr[face];
182  wTPtr[lPtr[sface]] -=
183  rDPtr[lPtr[sface]]*lowerPtr[sface]*wTPtr[uPtr[sface]];
184  }
185 }
186 
187 
188 // ************************************************************************* //
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
virtual void precondition(scalarField &wA, const scalarField &rA, const direction cmpt=0) const
Return wA the preconditioned form of residual rA.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
lduMatrix::preconditioner::addasymMatrixConstructorToTable< DILUPreconditioner > addDILUPreconditionerAsymMatrixConstructorToTable_
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
uint8_t direction
Definition: direction.H:45
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const lduMatrix & matrix() const
Definition: lduMatrix.H:226
scalarField & upper()
Definition: lduMatrix.C:197
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
DILUPreconditioner(const lduMatrix::solver &, const dictionary &solverControlsUnused)
Construct from matrix components and preconditioner solver controls.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:93
iterator begin()
Return an iterator to begin traversing the UList.
Definition: UListI.H:216
defineTypeNameAndDebug(combustionModel, 0)
virtual void preconditionT(scalarField &wT, const scalarField &rT, const direction cmpt=0) const
Return wT the transpose-matrix preconditioned form of residual rT.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
Abstract base-class for lduMatrix preconditioners.
Definition: lduMatrix.H:410
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:550
scalarField & lower()
Definition: lduMatrix.C:168
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
static void calcReciprocalD(scalarField &rD, const lduMatrix &matrix)
Calculate the reciprocal of the preconditioned diagonal.
scalarField & diag()
Definition: lduMatrix.C:186
Namespace for OpenFOAM.