PBiCGStab.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) 2016-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2021 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 "PBiCGStab.H"
30 #include "PrecisionAdaptor.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 
38  lduMatrix::solver::addsymMatrixConstructorToTable<PBiCGStab>
40 
41  lduMatrix::solver::addasymMatrixConstructorToTable<PBiCGStab>
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 Foam::PBiCGStab::PBiCGStab
49 (
50  const word& fieldName,
51  const lduMatrix& matrix,
52  const FieldField<Field, scalar>& interfaceBouCoeffs,
53  const FieldField<Field, scalar>& interfaceIntCoeffs,
54  const lduInterfaceFieldPtrsList& interfaces,
55  const dictionary& solverControls
56 )
57 :
59  (
60  fieldName,
61  matrix,
62  interfaceBouCoeffs,
63  interfaceIntCoeffs,
64  interfaces,
65  solverControls
66  )
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
73 (
75  const solveScalarField& source,
76  const direction cmpt
77 ) const
78 {
79  // --- Setup class containing solver performance data
80  solverPerformance solverPerf
81  (
82  lduMatrix::preconditioner::getName(controlDict_) + typeName,
83  fieldName_
84  );
85 
86  const label nCells = psi.size();
87 
88  solveScalar* __restrict__ psiPtr = psi.begin();
89 
90  solveScalarField pA(nCells);
91  solveScalar* __restrict__ pAPtr = pA.begin();
92 
93  solveScalarField yA(nCells);
94  solveScalar* __restrict__ yAPtr = yA.begin();
95 
96  // --- Calculate A.psi
97  matrix_.Amul(yA, psi, interfaceBouCoeffs_, interfaces_, cmpt);
98 
99  // --- Calculate initial residual field
100  solveScalarField rA(source - yA);
101  solveScalar* __restrict__ rAPtr = rA.begin();
102 
103  matrix().setResidualField
104  (
106  fieldName_,
107  true
108  );
109 
110  // --- Calculate normalisation factor
111  const solveScalar normFactor = this->normFactor(psi, source, yA, pA);
112 
113  if ((log_ >= 2) || (lduMatrix::debug >= 2))
114  {
115  Info<< " Normalisation factor = " << normFactor << endl;
116  }
117 
118  // --- Calculate normalised residual norm
119  solverPerf.initialResidual() =
120  gSumMag(rA, matrix().mesh().comm())
121  /normFactor;
122  solverPerf.finalResidual() = solverPerf.initialResidual();
123 
124  // --- Check convergence, solve if not converged
125  if
126  (
127  minIter_ > 0
128  || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
129  )
130  {
131  solveScalarField AyA(nCells);
132  solveScalar* __restrict__ AyAPtr = AyA.begin();
133 
134  solveScalarField sA(nCells);
135  solveScalar* __restrict__ sAPtr = sA.begin();
136 
137  solveScalarField zA(nCells);
138  solveScalar* __restrict__ zAPtr = zA.begin();
139 
140  solveScalarField tA(nCells);
141  solveScalar* __restrict__ tAPtr = tA.begin();
142 
143  // --- Store initial residual
144  const solveScalarField rA0(rA);
145 
146  // --- Initial values not used
147  solveScalar rA0rA = 0;
148  solveScalar alpha = 0;
149  solveScalar omega = 0;
150 
151  // --- Select and construct the preconditioner
152  if (!preconPtr_)
153  {
154  preconPtr_ = lduMatrix::preconditioner::New
155  (
156  *this,
157  controlDict_
158  );
159  }
160 
161  // --- Solver iteration
162  do
163  {
164  // --- Store previous rA0rA
165  const solveScalar rA0rAold = rA0rA;
166 
167  rA0rA = gSumProd(rA0, rA, matrix().mesh().comm());
168 
169  // --- Test for singularity
170  if (solverPerf.checkSingularity(mag(rA0rA)))
171  {
172  break;
173  }
174 
175  // --- Update pA
176  if (solverPerf.nIterations() == 0)
177  {
178  for (label cell=0; cell<nCells; cell++)
179  {
180  pAPtr[cell] = rAPtr[cell];
181  }
182  }
183  else
184  {
185  // --- Test for singularity
186  if (solverPerf.checkSingularity(mag(omega)))
187  {
188  break;
189  }
190 
191  const solveScalar beta = (rA0rA/rA0rAold)*(alpha/omega);
192 
193  for (label cell=0; cell<nCells; cell++)
194  {
195  pAPtr[cell] =
196  rAPtr[cell] + beta*(pAPtr[cell] - omega*AyAPtr[cell]);
197  }
198  }
199 
200  // --- Precondition pA
201  preconPtr_->precondition(yA, pA, cmpt);
202 
203  // --- Calculate AyA
204  matrix_.Amul(AyA, yA, interfaceBouCoeffs_, interfaces_, cmpt);
205 
206  const solveScalar rA0AyA =
207  gSumProd(rA0, AyA, matrix().mesh().comm());
208 
209  alpha = rA0rA/rA0AyA;
210 
211  // --- Calculate sA
212  for (label cell=0; cell<nCells; cell++)
213  {
214  sAPtr[cell] = rAPtr[cell] - alpha*AyAPtr[cell];
215  }
216 
217  // --- Test sA for convergence
218  solverPerf.finalResidual() =
219  gSumMag(sA, matrix().mesh().comm())/normFactor;
220 
221  if
222  (
223  solverPerf.nIterations() >= minIter_
224  && solverPerf.checkConvergence(tolerance_, relTol_, log_)
225  )
226  {
227  for (label cell=0; cell<nCells; cell++)
228  {
229  psiPtr[cell] += alpha*yAPtr[cell];
230  }
231 
232  solverPerf.nIterations()++;
233 
234  return solverPerf;
235  }
236 
237  // --- Precondition sA
238  preconPtr_->precondition(zA, sA, cmpt);
239 
240  // --- Calculate tA
241  matrix_.Amul(tA, zA, interfaceBouCoeffs_, interfaces_, cmpt);
242 
243  const solveScalar tAtA = gSumSqr(tA, matrix().mesh().comm());
244 
245  // --- Calculate omega from tA and sA
246  // (cheaper than using zA with preconditioned tA)
247  omega = gSumProd(tA, sA, matrix().mesh().comm())/tAtA;
248 
249  // --- Update solution and residual
250  for (label cell=0; cell<nCells; cell++)
251  {
252  psiPtr[cell] += alpha*yAPtr[cell] + omega*zAPtr[cell];
253  rAPtr[cell] = sAPtr[cell] - omega*tAPtr[cell];
254  }
255 
256  solverPerf.finalResidual() =
257  gSumMag(rA, matrix().mesh().comm())
258  /normFactor;
259  } while
260  (
261  (
262  ++solverPerf.nIterations() < maxIter_
263  && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
264  )
265  || solverPerf.nIterations() < minIter_
266  );
267  }
268 
269  if (preconPtr_)
270  {
271  preconPtr_->setFinished(solverPerf);
272  }
273 
274  matrix().setResidualField
275  (
276  ConstPrecisionAdaptor<scalar, solveScalar>(rA)(),
277  fieldName_,
278  false
279  );
280 
281  return solverPerf;
282 }
283 
284 
286 (
287  scalarField& psi_s,
288  const scalarField& source,
289  const direction cmpt
290 ) const
291 {
293  return scalarSolve
294  (
295  tpsi.ref(),
297  cmpt
298  );
299 }
300 
301 
302 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:46
Base solver class.
Definition: solver.H:45
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< solveScalar > solveScalarField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void setResidualField(const scalarField &residual, const word &fieldName, const bool initial) const
Set the residual field using an IOField on the object registry if it exists.
Definition: lduMatrix.C:327
scalarProduct< Type, Type >::type gSumProd(const UList< Type > &f1, const UList< Type > &f2, const label comm)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A non-const Field/List wrapper with possible data conversion.
lduMatrix::solver::addsymMatrixConstructorToTable< PBiCGStab > addPBiCGStabSymMatrixConstructorToTable_
Definition: PBiCGStab.C:32
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
static autoPtr< preconditioner > New(const solver &sol, const dictionary &solverControls)
Return a new preconditioner.
A const Field/List wrapper with possible data conversion.
lduMatrix::solver::addasymMatrixConstructorToTable< PBiCGStab > addPBiCGStabAsymMatrixConstructorToTable_
Definition: PBiCGStab.C:35
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:391
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition: PBiCGStab.C:279
int debug
Static debugging option.
Container< Type > & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: refPtrI.H:230
defineTypeNameAndDebug(combustionModel, 0)
outerProduct1< Type >::type gSumSqr(const UList< Type > &f, const label comm)
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
static word getName(const dictionary &)
Find the preconditioner name (directly or from a sub-dictionary)
const volScalarField & psi
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Preconditioned bi-conjugate gradient stabilized solver for asymmetric lduMatrices using a run-time se...
Definition: PBiCGStab.H:63
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition: PBiCGStab.C:66
Namespace for OpenFOAM.