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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2023 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 "PCG.H"
30 #include "PrecisionAdaptor.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 
38  lduMatrix::solver::addsymMatrixConstructorToTable<PCG>
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::PCG::PCG
46 (
47  const word& fieldName,
48  const lduMatrix& matrix,
49  const FieldField<Field, scalar>& interfaceBouCoeffs,
50  const FieldField<Field, scalar>& interfaceIntCoeffs,
51  const lduInterfaceFieldPtrsList& interfaces,
52  const dictionary& solverControls
53 )
54 :
56  (
57  fieldName,
58  matrix,
59  interfaceBouCoeffs,
60  interfaceIntCoeffs,
61  interfaces,
62  solverControls
63  )
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
70 (
72  const solveScalarField& source,
73  const direction cmpt
74 ) const
75 {
76  // --- Setup class containing solver performance data
77  solverPerformance solverPerf
78  (
79  lduMatrix::preconditioner::getName(controlDict_) + typeName,
80  fieldName_
81  );
82 
83  label nCells = psi.size();
84 
85  solveScalar* __restrict__ psiPtr = psi.begin();
86 
87  solveScalarField pA(nCells);
88  solveScalar* __restrict__ pAPtr = pA.begin();
89 
90  solveScalarField wA(nCells);
91  solveScalar* __restrict__ wAPtr = wA.begin();
92 
93  solveScalar wArA = solverPerf.great_;
94  solveScalar wArAold = wArA;
95 
96  // --- Calculate A.psi
97  matrix_.Amul(wA, psi, interfaceBouCoeffs_, interfaces_, cmpt);
98 
99  // --- Calculate initial residual field
100  solveScalarField rA(source - wA);
101  solveScalar* __restrict__ rAPtr = rA.begin();
102 
103  matrix().setResidualField
104  (
106  fieldName_,
107  true
108  );
109 
110  // --- Calculate normalisation factor
111  solveScalar normFactor = this->normFactor(psi, source, wA, 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  // --- Select and construct the preconditioner
132  if (!preconPtr_)
133  {
134  preconPtr_ = lduMatrix::preconditioner::New
135  (
136  *this,
137  controlDict_
138  );
139  }
140 
141  // --- Solver iteration
142  do
143  {
144  // --- Store previous wArA
145  wArAold = wArA;
146 
147  // --- Precondition residual
148  preconPtr_->precondition(wA, rA, cmpt);
149 
150  // --- Update search directions:
151  wArA = gSumProd(wA, rA, matrix().mesh().comm());
152 
153  if (solverPerf.nIterations() == 0)
154  {
155  for (label cell=0; cell<nCells; cell++)
156  {
157  pAPtr[cell] = wAPtr[cell];
158  }
159  }
160  else
161  {
162  const solveScalar beta = wArA/wArAold;
163 
164  for (label cell=0; cell<nCells; cell++)
165  {
166  pAPtr[cell] = wAPtr[cell] + beta*pAPtr[cell];
167  }
168  }
169 
170 
171  // --- Update preconditioned residual
172  matrix_.Amul(wA, pA, interfaceBouCoeffs_, interfaces_, cmpt);
173 
174  solveScalar wApA = gSumProd(wA, pA, matrix().mesh().comm());
175 
176  // --- Test for singularity
177  if (solverPerf.checkSingularity(mag(wApA)/normFactor)) break;
178 
179 
180  // --- Update solution and residual:
181 
182  const solveScalar alpha = wArA/wApA;
183 
184  for (label cell=0; cell<nCells; cell++)
185  {
186  psiPtr[cell] += alpha*pAPtr[cell];
187  rAPtr[cell] -= alpha*wAPtr[cell];
188  }
189 
190  solverPerf.finalResidual() =
191  gSumMag(rA, matrix().mesh().comm())
192  /normFactor;
193 
194  } while
195  (
196  (
197  ++solverPerf.nIterations() < maxIter_
198  && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
199  )
200  || solverPerf.nIterations() < minIter_
201  );
202  }
203 
204  if (preconPtr_)
205  {
206  preconPtr_->setFinished(solverPerf);
207  }
208 
209  matrix().setResidualField
210  (
211  ConstPrecisionAdaptor<scalar, solveScalar>(rA)(),
212  fieldName_,
213  false
214  );
215 
216  return solverPerf;
217 }
218 
219 
221 (
222  scalarField& psi_s,
223  const scalarField& source,
224  const direction cmpt
225 ) const
226 {
228  return scalarSolve
229  (
230  tpsi.ref(),
232  cmpt
233  );
234 }
235 
236 
237 // ************************************************************************* //
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)
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)
lduMatrix::solver::addsymMatrixConstructorToTable< PCG > addPCGSymMatrixConstructorToTable_
Definition: PCG.C:32
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A non-const Field/List wrapper with possible data conversion.
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition: PCG.C:214
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.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:391
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)
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition: PCG.C:63
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 conjugate gradient solver for symmetric lduMatrices using a run-time selectable precon...
Definition: PCG.H:51
Namespace for OpenFOAM.