PPCG.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) 2019-2020 Mattijs Janssens
9  Copyright (C) 2020-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 "PPCG.H"
30 #include "PrecisionAdaptor.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 
38  lduMatrix::solver::addsymMatrixConstructorToTable<PPCG>
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::PPCG::gSumMagProd
46 (
47  FixedList<solveScalar, 3>& globalSum,
48  const solveScalarField& a,
49  const solveScalarField& b,
50  const solveScalarField& c,
51  const solveScalarField& sumMag,
52  UPstream::Request& request,
53  const label comm
54 )
55 {
56  const label nCells = a.size();
57 
58  globalSum = 0.0;
59  for (label cell=0; cell<nCells; ++cell)
60  {
61  globalSum[0] += a[cell]*b[cell]; // sumProd(a, b)
62  globalSum[1] += a[cell]*c[cell]; // sumProd(a, c)
63  globalSum[2] += mag(sumMag[cell]);
64  }
65 
66  if (UPstream::parRun())
67  {
69  (
70  globalSum.data(),
71  globalSum.size(),
72  sumOp<solveScalar>(),
73  UPstream::msgType(), // (ignored): direct MPI call
74  comm,
75  request
76  );
77  }
78 }
79 
80 
82 (
84  const solveScalarField& source,
85  const direction cmpt,
86  const bool cgMode
87 ) const
88 {
89  // --- Setup class containing solver performance data
90  solverPerformance solverPerf
91  (
92  lduMatrix::preconditioner::getName(controlDict_) + type(),
93  fieldName_
94  );
95 
96  const label comm = matrix().mesh().comm();
97  const label nCells = psi.size();
98  solveScalarField w(nCells);
99 
100  // --- Calculate A.psi
101  matrix_.Amul(w, psi, interfaceBouCoeffs_, interfaces_, cmpt);
102 
103  // --- Calculate initial residual field
104  solveScalarField r(source - w);
105 
106  // --- Calculate normalisation factor
107  solveScalarField p(nCells);
108  const solveScalar normFactor = this->normFactor(psi, source, w, p);
109 
110  if ((log_ >= 2) || (lduMatrix::debug >= 2))
111  {
112  Info<< " Normalisation factor = " << normFactor << endl;
113  }
114 
115  // --- Select and construct the preconditioner
116  if (!preconPtr_)
117  {
118  preconPtr_ = lduMatrix::preconditioner::New
119  (
120  *this,
121  controlDict_
122  );
123  }
124 
125  // --- Precondition residual (= u0)
126  solveScalarField u(nCells);
127  preconPtr_->precondition(u, r, cmpt);
128 
129  // --- Calculate A*u - reuse w
130  matrix_.Amul(w, u, interfaceBouCoeffs_, interfaces_, cmpt);
131 
132 
133  // State
134  solveScalarField s(nCells);
135  solveScalarField q(nCells);
136  solveScalarField z(nCells);
137 
138  solveScalarField m(nCells);
139 
140  FixedList<solveScalar, 3> globalSum;
141  UPstream::Request outstandingRequest;
142  if (cgMode)
143  {
144  // --- Start global reductions for inner products
145  gSumMagProd(globalSum, u, r, w, r, outstandingRequest, comm);
146 
147  // --- Precondition residual
148  preconPtr_->precondition(m, w, cmpt);
149  }
150  else
151  {
152  // --- Precondition residual
153  preconPtr_->precondition(m, w, cmpt);
154 
155  // --- Start global reductions for inner products
156  gSumMagProd(globalSum, w, u, m, r, outstandingRequest, comm);
157  }
158 
159  // --- Calculate A*m
160  solveScalarField n(nCells);
161  matrix_.Amul(n, m, interfaceBouCoeffs_, interfaces_, cmpt);
162 
163  solveScalar alpha = 0.0;
164  solveScalar gamma = 0.0;
165 
166  // --- Solver iteration
167  for
168  (
169  solverPerf.nIterations() = 0;
170  solverPerf.nIterations() < maxIter_;
171  solverPerf.nIterations()++
172  )
173  {
174  // Make sure gamma,delta are available
175  outstandingRequest.wait();
176 
177  const solveScalar gammaOld = gamma;
178  gamma = globalSum[0];
179  const solveScalar delta = globalSum[1];
180 
181  solverPerf.finalResidual() = globalSum[2]/normFactor;
182  if (solverPerf.nIterations() == 0)
183  {
184  solverPerf.initialResidual() = solverPerf.finalResidual();
185  }
186 
187  // Check convergence (bypass if not enough iterations yet)
188  if
189  (
190  (minIter_ <= 0 || solverPerf.nIterations() >= minIter_)
191  && solverPerf.checkConvergence(tolerance_, relTol_, log_)
192  )
193  {
194  break;
195  }
196 
197 
198  if (solverPerf.nIterations() == 0)
199  {
200  alpha = gamma/delta;
201  z = n;
202  q = m;
203  s = w;
204  p = u;
205  }
206  else
207  {
208  const solveScalar beta = gamma/gammaOld;
210 
211  for (label cell=0; cell<nCells; ++cell)
212  {
213  z[cell] = n[cell] + beta*z[cell];
214  q[cell] = m[cell] + beta*q[cell];
215  s[cell] = w[cell] + beta*s[cell];
216  p[cell] = u[cell] + beta*p[cell];
217  }
218  }
219 
220  for (label cell=0; cell<nCells; ++cell)
221  {
222  psi[cell] += alpha*p[cell];
223  r[cell] -= alpha*s[cell];
224  u[cell] -= alpha*q[cell];
225  w[cell] -= alpha*z[cell];
226  }
227 
228  if (cgMode)
229  {
230  // --- Start global reductions for inner products
231  gSumMagProd(globalSum, u, r, w, r, outstandingRequest, comm);
232 
233  // --- Precondition residual
234  preconPtr_->precondition(m, w, cmpt);
235  }
236  else
237  {
238  // --- Precondition residual
239  preconPtr_->precondition(m, w, cmpt);
240 
241  // --- Start global reductions for inner products
242  gSumMagProd(globalSum, w, u, m, r, outstandingRequest, comm);
243  }
244 
245  // --- Calculate A*m
246  matrix_.Amul(n, m, interfaceBouCoeffs_, interfaces_, cmpt);
247  }
248 
249  // Cleanup any outstanding requests
250  outstandingRequest.wait();
251 
252  if (preconPtr_)
253  {
254  preconPtr_->setFinished(solverPerf);
255  }
256 
257  //TBD
258  //matrix().setResidualField
259  //(
260  // ConstPrecisionAdaptor<scalar, solveScalar>(rA)(),
261  // fieldName_,
262  // false
263  //);
264 
265  return solverPerf;
266 }
267 
268 
269 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
270 
271 Foam::PPCG::PPCG
272 (
273  const word& fieldName,
274  const lduMatrix& matrix,
275  const FieldField<Field, scalar>& interfaceBouCoeffs,
276  const FieldField<Field, scalar>& interfaceIntCoeffs,
277  const lduInterfaceFieldPtrsList& interfaces,
278  const dictionary& solverControls
279 )
280 :
282  (
283  fieldName,
284  matrix,
285  interfaceBouCoeffs,
286  interfaceIntCoeffs,
287  interfaces,
288  solverControls
289  )
290 {}
291 
292 
293 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
294 
296 (
297  scalarField& psi_s,
298  const scalarField& source,
299  const direction cmpt
300 ) const
301 {
303  return scalarSolveCG
304  (
305  tpsi.ref(),
307  cmpt,
308  true // operate in conjugate-gradient mode
309  );
310 }
311 
312 
313 // ************************************************************************* //
scalar delta
dimensioned< typename typeOfMag< Type >::type > sumMag(const DimensionedField< Type, GeoMesh > &f1)
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
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve the matrix with this solver.
Definition: PPCG.C:289
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
virtual label comm() const =0
Return communicator used for parallel communication.
A non-const Field/List wrapper with possible data conversion.
Preconditioned pipelined conjugate gradient solver for symmetric lduMatrices using a run-time selecta...
Definition: PPCG.H:60
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
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.
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition: lduMatrix.H:718
const dimensionedScalar c
Speed of light in a vacuum.
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)
label n
static word getName(const dictionary &)
Find the preconditioner name (directly or from a sub-dictionary)
const volScalarField & psi
const scalar gamma
Definition: EEqn.H:9
volScalarField & p
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
lduMatrix::solver::addsymMatrixConstructorToTable< PPCG > addPPCGSymMatrixConstructorToTable_
Definition: PPCG.C:32
Namespace for OpenFOAM.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
solverPerformance scalarSolveCG(solveScalarField &psi, const solveScalarField &source, const direction cmpt, const bool cgMode) const
CG solver. Operates either in conjugate-gradient mode or conjugate residual.
Definition: PPCG.C:75