FreeFOAM The Cross-Platform CFD Toolkit
PCG.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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 <OpenFOAM/PCG.H>
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(PCG, 0);
33 
34  lduMatrix::solver::addsymMatrixConstructorToTable<PCG>
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 Foam::PCG::PCG
42 (
43  const word& fieldName,
44  const lduMatrix& matrix,
45  const FieldField<Field, scalar>& interfaceBouCoeffs,
46  const FieldField<Field, scalar>& interfaceIntCoeffs,
47  const lduInterfaceFieldPtrsList& interfaces,
48  const dictionary& solverControls
49 )
50 :
52  (
53  fieldName,
54  matrix,
55  interfaceBouCoeffs,
56  interfaceIntCoeffs,
57  interfaces,
58  solverControls
59  )
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
66 (
68  const scalarField& source,
69  const direction cmpt
70 ) const
71 {
72  // --- Setup class containing solver performance data
74  (
75  lduMatrix::preconditioner::getName(controlDict_) + typeName,
76  fieldName_
77  );
78 
79  register label nCells = psi.size();
80 
81  scalar* __restrict__ psiPtr = psi.begin();
82 
83  scalarField pA(nCells);
84  scalar* __restrict__ pAPtr = pA.begin();
85 
86  scalarField wA(nCells);
87  scalar* __restrict__ wAPtr = wA.begin();
88 
89  scalar wArA = matrix_.great_;
90  scalar wArAold = wArA;
91 
92  // --- Calculate A.psi
93  matrix_.Amul(wA, psi, interfaceBouCoeffs_, interfaces_, cmpt);
94 
95  // --- Calculate initial residual field
96  scalarField rA(source - wA);
97  scalar* __restrict__ rAPtr = rA.begin();
98 
99  // --- Calculate normalisation factor
100  scalar normFactor = this->normFactor(psi, source, wA, pA);
101 
102  if (lduMatrix::debug >= 2)
103  {
104  Info<< " Normalisation factor = " << normFactor << endl;
105  }
106 
107  // --- Calculate normalised residual norm
108  solverPerf.initialResidual() = gSumMag(rA)/normFactor;
109  solverPerf.finalResidual() = solverPerf.initialResidual();
110 
111  // --- Check convergence, solve if not converged
112  if (!solverPerf.checkConvergence(tolerance_, relTol_))
113  {
114  // --- Select and construct the preconditioner
116  lduMatrix::preconditioner::New
117  (
118  *this,
119  controlDict_
120  );
121 
122  // --- Solver iteration
123  do
124  {
125  // --- Store previous wArA
126  wArAold = wArA;
127 
128  // --- Precondition residual
129  preconPtr->precondition(wA, rA, cmpt);
130 
131  // --- Update search directions:
132  wArA = gSumProd(wA, rA);
133 
134  if (solverPerf.nIterations() == 0)
135  {
136  for (register label cell=0; cell<nCells; cell++)
137  {
138  pAPtr[cell] = wAPtr[cell];
139  }
140  }
141  else
142  {
143  scalar beta = wArA/wArAold;
144 
145  for (register label cell=0; cell<nCells; cell++)
146  {
147  pAPtr[cell] = wAPtr[cell] + beta*pAPtr[cell];
148  }
149  }
150 
151 
152  // --- Update preconditioned residual
153  matrix_.Amul(wA, pA, interfaceBouCoeffs_, interfaces_, cmpt);
154 
155  scalar wApA = gSumProd(wA, pA);
156 
157 
158  // --- Test for singularity
159  if (solverPerf.checkSingularity(mag(wApA)/normFactor)) break;
160 
161 
162  // --- Update solution and residual:
163 
164  scalar alpha = wArA/wApA;
165 
166  for (register label cell=0; cell<nCells; cell++)
167  {
168  psiPtr[cell] += alpha*pAPtr[cell];
169  rAPtr[cell] -= alpha*wAPtr[cell];
170  }
171 
172  solverPerf.finalResidual() = gSumMag(rA)/normFactor;
173 
174  } while
175  (
176  solverPerf.nIterations()++ < maxIter_
177  && !(solverPerf.checkConvergence(tolerance_, relTol_))
178  );
179  }
180 
181  return solverPerf;
182 }
183 
184 
185 // ************************ vim: set sw=4 sts=4 et: ************************ //