DILUPreconditioner.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-2015 OpenFOAM Foundation
9  Copyright (C) 2019 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 "DILUPreconditioner.H"
30 #include <algorithm>
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(DILUPreconditioner, 0);
37 
38  lduMatrix::preconditioner::
39  addasymMatrixConstructorToTable<DILUPreconditioner>
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const lduMatrix::solver& sol,
49  const dictionary&
50 )
51 :
53  rD_(sol.matrix().diag().size())
54 {
55  const scalarField& diag = sol.matrix().diag();
56  std::copy(diag.begin(), diag.end(), rD_.begin());
57 
58  calcReciprocalD(rD_, sol.matrix());
59 }
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
65 (
66  solveScalarField& rD,
67  const lduMatrix& matrix
68 )
69 {
70  solveScalar* __restrict__ rDPtr = rD.begin();
71 
72  const label* const __restrict__ uPtr = matrix.lduAddr().upperAddr().begin();
73  const label* const __restrict__ lPtr = matrix.lduAddr().lowerAddr().begin();
74 
75  const scalar* const __restrict__ upperPtr = matrix.upper().begin();
76  const scalar* const __restrict__ lowerPtr = matrix.lower().begin();
77 
78  label nFaces = matrix.upper().size();
79  for (label face=0; face<nFaces; face++)
80  {
81  rDPtr[uPtr[face]] -= upperPtr[face]*lowerPtr[face]/rDPtr[lPtr[face]];
82  }
83 
84 
85  // Calculate the reciprocal of the preconditioned diagonal
86  const label nCells = rD.size();
87 
88  for (label cell=0; cell<nCells; cell++)
89  {
90  rDPtr[cell] = 1.0/rDPtr[cell];
91  }
92 }
93 
94 
96 (
97  solveScalarField& wA,
98  const solveScalarField& rA,
99  const direction
100 ) const
101 {
102  solveScalar* __restrict__ wAPtr = wA.begin();
103  const solveScalar* __restrict__ rAPtr = rA.begin();
104  const solveScalar* __restrict__ rDPtr = rD_.begin();
105 
106  const label* const __restrict__ uPtr =
107  solver_.matrix().lduAddr().upperAddr().begin();
108  const label* const __restrict__ lPtr =
109  solver_.matrix().lduAddr().lowerAddr().begin();
110  const label* const __restrict__ losortPtr =
111  solver_.matrix().lduAddr().losortAddr().begin();
112 
113  const scalar* const __restrict__ upperPtr =
114  solver_.matrix().upper().begin();
115  const scalar* const __restrict__ lowerPtr =
116  solver_.matrix().lower().begin();
117 
118  const label nCells = wA.size();
119  const label nFaces = solver_.matrix().upper().size();
120  const label nFacesM1 = nFaces - 1;
121 
122  for (label cell=0; cell<nCells; cell++)
123  {
124  wAPtr[cell] = rDPtr[cell]*rAPtr[cell];
125  }
126 
127  for (label face=0; face<nFaces; face++)
128  {
129  const label sface = losortPtr[face];
130  wAPtr[uPtr[sface]] -=
131  rDPtr[uPtr[sface]]*lowerPtr[sface]*wAPtr[lPtr[sface]];
132  }
133 
134  for (label face=nFacesM1; face>=0; face--)
135  {
136  wAPtr[lPtr[face]] -=
137  rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[uPtr[face]];
138  }
139 }
140 
141 
143 (
144  solveScalarField& wT,
145  const solveScalarField& rT,
146  const direction
147 ) const
148 {
149  solveScalar* __restrict__ wTPtr = wT.begin();
150  const solveScalar* __restrict__ rTPtr = rT.begin();
151  const solveScalar* __restrict__ rDPtr = rD_.begin();
152 
153  const label* const __restrict__ uPtr =
154  solver_.matrix().lduAddr().upperAddr().begin();
155  const label* const __restrict__ lPtr =
156  solver_.matrix().lduAddr().lowerAddr().begin();
157  const label* const __restrict__ losortPtr =
158  solver_.matrix().lduAddr().losortAddr().begin();
159 
160  const scalar* const __restrict__ upperPtr =
161  solver_.matrix().upper().begin();
162  const scalar* const __restrict__ lowerPtr =
163  solver_.matrix().lower().begin();
164 
165  const label nCells = wT.size();
166  const label nFaces = solver_.matrix().upper().size();
167  const label nFacesM1 = nFaces - 1;
168 
169  for (label cell=0; cell<nCells; cell++)
170  {
171  wTPtr[cell] = rDPtr[cell]*rTPtr[cell];
172  }
173 
174  for (label face=0; face<nFaces; face++)
175  {
176  wTPtr[uPtr[face]] -=
177  rDPtr[uPtr[face]]*upperPtr[face]*wTPtr[lPtr[face]];
178  }
179 
180 
181  for (label face=nFacesM1; face>=0; face--)
182  {
183  const label sface = losortPtr[face];
184  wTPtr[lPtr[sface]] -=
185  rDPtr[lPtr[sface]]*lowerPtr[sface]*wTPtr[uPtr[sface]];
186  }
187 }
188 
189 
190 // ************************************************************************* //
Foam::lduMatrix::lduAddr
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:578
Foam::DILUPreconditioner::preconditionT
virtual void preconditionT(solveScalarField &wT, const solveScalarField &rT, const direction cmpt=0) const
Return wT the transpose-matrix preconditioned form of residual rT.
Definition: DILUPreconditioner.C:143
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:287
Foam::lduMatrix::upper
scalarField & upper()
Definition: lduMatrix.C:203
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:83
Foam::lduMatrix::solver
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:98
Foam::DILUPreconditioner::DILUPreconditioner
DILUPreconditioner(const lduMatrix::solver &, const dictionary &solverControlsUnused)
Construct from matrix components and preconditioner solver controls.
Definition: DILUPreconditioner.C:47
Foam::lduAddressing::upperAddr
virtual const labelUList & upperAddr() const =0
Return upper addressing.
Foam::DILUPreconditioner::precondition
virtual void precondition(solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const
Return wA the preconditioned form of residual rA.
Definition: DILUPreconditioner.C:96
Foam::UList::begin
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:329
Foam::Field< scalar >
Foam::DILUPreconditioner::calcReciprocalD
static void calcReciprocalD(solveScalarField &, const lduMatrix &)
Calculate the reciprocal of the preconditioned diagonal.
Definition: DILUPreconditioner.C:65
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::lduMatrix::diag
scalarField & diag()
Definition: lduMatrix.C:192
Foam::lduAddressing::lowerAddr
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Foam::addDILUPreconditionerAsymMatrixConstructorToTable_
lduMatrix::preconditioner::addasymMatrixConstructorToTable< DILUPreconditioner > addDILUPreconditionerAsymMatrixConstructorToTable_
Definition: DILUPreconditioner.C:40
Foam::direction
uint8_t direction
Definition: direction.H:52
DILUPreconditioner.H
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::lduMatrix::preconditioner
Abstract base-class for lduMatrix preconditioners.
Definition: lduMatrix.H:433
Foam::lduMatrix::lower
scalarField & lower()
Definition: lduMatrix.C:174
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::lduMatrix::solver::matrix
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:233