FPCG.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-2021 OpenCFD Ltd.
10  Copyright (C) 2023 Huawei (Yu Ankun)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "FPCG.H"
31 #include "PrecisionAdaptor.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38 
39  lduMatrix::solver::addsymMatrixConstructorToTable<FPCG>
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 Foam::FPCG::FPCG
47 (
48  const word& fieldName,
49  const lduMatrix& matrix,
50  const FieldField<Field, scalar>& interfaceBouCoeffs,
51  const FieldField<Field, scalar>& interfaceIntCoeffs,
52  const lduInterfaceFieldPtrsList& interfaces,
53  const dictionary& solverControls
54 )
55 :
57  (
58  fieldName,
59  matrix,
60  interfaceBouCoeffs,
61  interfaceIntCoeffs,
62  interfaces,
63  solverControls
64  )
65 {}
66 
67 
68 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
69 
70 void Foam::FPCG::gSumMagProd
71 (
72  FixedList<solveScalar, 2>& globalSum,
73  const solveScalarField& a,
74  const solveScalarField& b,
75  const label comm
76 )
77 {
78  const label nCells = a.size();
79 
80  globalSum = 0.0;
81  for (label cell=0; cell<nCells; ++cell)
82  {
83  globalSum[0] += a[cell]*b[cell]; // sumProd(a, b)
84  globalSum[1] += mag(b[cell]); // sumMag(b)
85  }
86 
88  (
89  globalSum.data(),
90  globalSum.size(),
91  sumOp<solveScalar>(),
92  UPstream::msgType(), // (ignored): direct MPI call
93  comm
94  );
95 }
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
100 (
102  const solveScalarField& source,
103  const direction cmpt
104 ) const
105 {
106  // --- Setup class containing solver performance data
107  solverPerformance solverPerf
108  (
109  lduMatrix::preconditioner::getName(controlDict_) + typeName,
110  fieldName_
111  );
112 
113  label nCells = psi.size();
114 
115  solveScalar* __restrict__ psiPtr = psi.begin();
116 
117  solveScalarField pA(nCells);
118  solveScalar* __restrict__ pAPtr = pA.begin();
119 
120  solveScalarField wA(nCells);
121  solveScalar* __restrict__ wAPtr = wA.begin();
122 
123  solveScalar wArA = solverPerf.great_;
124  solveScalar wArAold = wArA;
125 
126  // --- Calculate A.psi
127  matrix_.Amul(wA, psi, interfaceBouCoeffs_, interfaces_, cmpt);
128 
129  // --- Calculate initial residual field
130  solveScalarField rA(source - wA);
131  solveScalar* __restrict__ rAPtr = rA.begin();
132 
133  matrix().setResidualField
134  (
136  fieldName_,
137  true
138  );
139 
140  // --- Calculate normalisation factor
141  solveScalar normFactor = this->normFactor(psi, source, wA, pA);
142 
143  if ((log_ >= 2) || (lduMatrix::debug >= 2))
144  {
145  Info<< " Normalisation factor = " << normFactor << endl;
146  }
147 
148  // --- Calculate normalised residual norm
149  solverPerf.initialResidual() =
150  gSumMag(rA, matrix().mesh().comm())
151  /normFactor;
152  solverPerf.finalResidual() = solverPerf.initialResidual();
153 
154  // --- Check convergence, solve if not converged
155  if
156  (
157  minIter_ > 0
158  || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
159  )
160  {
161  // --- Select and construct the preconditioner
162  if (!preconPtr_)
163  {
164  preconPtr_ = lduMatrix::preconditioner::New
165  (
166  *this,
167  controlDict_
168  );
169  }
170 
171  FixedList<solveScalar, 2> globalSum;
172 
173  // --- Solver iteration
174  do
175  {
176  // --- Store previous wArA
177  wArAold = wArA;
178 
179  // --- Precondition residual
180  preconPtr_->precondition(wA, rA, cmpt);
181 
182  // --- Update search directions and calculate residual:
183  gSumMagProd(globalSum, wA, rA, matrix().mesh().comm());
184 
185  wArA = globalSum[0];
186 
187  solverPerf.finalResidual() = globalSum[1]/normFactor;
188 
189  // Check convergence (bypass if not enough iterations yet)
190  if
191  (
192  (minIter_ <= 0 || solverPerf.nIterations() >= minIter_)
193  && solverPerf.checkConvergence(tolerance_, relTol_, log_)
194  )
195  {
196  break;
197  }
198 
199  if (solverPerf.nIterations() == 0)
200  {
201  for (label cell=0; cell<nCells; cell++)
202  {
203  pAPtr[cell] = wAPtr[cell];
204  }
205  }
206  else
207  {
208  const solveScalar beta = wArA/wArAold;
209 
210  for (label cell=0; cell<nCells; cell++)
211  {
212  pAPtr[cell] = wAPtr[cell] + beta*pAPtr[cell];
213  }
214  }
215 
216 
217  // --- Update preconditioned residual
218  matrix_.Amul(wA, pA, interfaceBouCoeffs_, interfaces_, cmpt);
219 
220  solveScalar wApA = gSumProd(wA, pA, matrix().mesh().comm());
221 
222  // --- Test for singularity
223  if (solverPerf.checkSingularity(mag(wApA)/normFactor)) break;
224 
225 
226  // --- Update solution and residual:
227 
228  const solveScalar alpha = wArA/wApA;
229 
230  for (label cell=0; cell<nCells; cell++)
231  {
232  psiPtr[cell] += alpha*pAPtr[cell];
233  rAPtr[cell] -= alpha*wAPtr[cell];
234  }
235 
236  } while
237  (
238  ++solverPerf.nIterations() < maxIter_
239  );
240  }
241 
242  if (preconPtr_)
243  {
244  preconPtr_->setFinished(solverPerf);
245  }
246 
247  matrix().setResidualField
248  (
249  ConstPrecisionAdaptor<scalar, solveScalar>(rA)(),
250  fieldName_,
251  false
252  );
254  return solverPerf;
255 }
256 
257 
258 
260 (
261  scalarField& psi_s,
262  const scalarField& source,
263  const direction cmpt
264 ) const
265 {
267  return scalarSolve
268  (
269  tpsi.ref(),
271  cmpt
272  );
273 }
274 
275 
276 // ************************************************************************* //
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)
virtual solverPerformance scalarSolve(solveScalarField &psi, const solveScalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition: FPCG.C:93
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
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition: FPCG.C:253
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
A non-const Field/List wrapper with possible data conversion.
bool checkSingularity(const Type &residual)
Singularity test.
bool checkConvergence(const Type &tolerance, const Type &relTolerance, const int logLevel=0)
Check, store and return convergence.
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
static const scalar great_
Large Type for the use in solvers.
dynamicFvMesh & mesh
const labelType & nIterations() const noexcept
Return number of iterations.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
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.
const Type & finalResidual() const noexcept
Return final residual.
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)
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
lduMatrix::solver::addsymMatrixConstructorToTable< FPCG > addFPCGSymMatrixConstructorToTable_
Definition: FPCG.C:33
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
A &#39;faster&#39; preconditioned conjugate gradient solver for symmetric lduMatrices using a run-time select...
Definition: FPCG.H:55
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: [].
Namespace for OpenFOAM.
const Type & initialResidual() const noexcept
Return initial residual.