DICPreconditioner.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 "DICPreconditioner.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(DICPreconditioner, 0);
33 
34  lduMatrix::preconditioner::
35  addsymMatrixConstructorToTable<DICPreconditioner>
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  const scalar* const __restrict__ upperPtr = matrix.upper().begin();
68 
69  // Calculate the DIC diagonal
70  const label nFaces = matrix.upper().size();
71  for (label face=0; face<nFaces; face++)
72  {
73  rDPtr[uPtr[face]] -= upperPtr[face]*upperPtr[face]/rDPtr[lPtr[face]];
74  }
75 
76 
77  // Calculate the reciprocal of the preconditioned diagonal
78  const label nCells = rD.size();
79 
80  for (label cell=0; cell<nCells; cell++)
81  {
82  rDPtr[cell] = 1.0/rDPtr[cell];
83  }
84 }
85 
86 
88 (
89  scalarField& wA,
90  const scalarField& rA,
91  const direction
92 ) const
93 {
94  scalar* __restrict__ wAPtr = wA.begin();
95  const scalar* __restrict__ rAPtr = rA.begin();
96  const scalar* __restrict__ rDPtr = rD_.begin();
97 
98  const label* const __restrict__ uPtr =
99  solver_.matrix().lduAddr().upperAddr().begin();
100  const label* const __restrict__ lPtr =
101  solver_.matrix().lduAddr().lowerAddr().begin();
102  const scalar* const __restrict__ upperPtr =
103  solver_.matrix().upper().begin();
104 
105  label nCells = wA.size();
106  label nFaces = solver_.matrix().upper().size();
107  label nFacesM1 = nFaces - 1;
108 
109  for (label cell=0; cell<nCells; cell++)
110  {
111  wAPtr[cell] = rDPtr[cell]*rAPtr[cell];
112  }
113 
114  for (label face=0; face<nFaces; face++)
115  {
116  wAPtr[uPtr[face]] -= rDPtr[uPtr[face]]*upperPtr[face]*wAPtr[lPtr[face]];
117  }
118 
119  for (label face=nFacesM1; face>=0; face--)
120  {
121  wAPtr[lPtr[face]] -= rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[uPtr[face]];
122  }
123 }
124 
125 
126 // ************************************************************************* //
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
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.
lduMatrix::preconditioner::addsymMatrixConstructorToTable< DICPreconditioner > addDICPreconditionerSymMatrixConstructorToTable_
DICPreconditioner(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
virtual void precondition(scalarField &wA, const scalarField &rA, const direction cmpt=0) const
Return wA the preconditioned form of residual rA.
defineTypeNameAndDebug(combustionModel, 0)
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
static void calcReciprocalD(scalarField &rD, const lduMatrix &matrix)
Calculate the reciprocal of the preconditioned diagonal.
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:550
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
scalarField & diag()
Definition: lduMatrix.C:186
Namespace for OpenFOAM.