GAMGSolverSolve.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) 2016-2021,2023 OpenCFD Ltd.
10  Copyright (C) 2023 Huawei (Yu Ankun)
11  Copyright (C) 2023 OpenCFD Ltd.
12 -------------------------------------------------------------------------------
13 License
14  This file is part of OpenFOAM.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "GAMGSolver.H"
32 #include "SubField.H"
33 #include "PrecisionAdaptor.H"
34 
35 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36 
38 (
39  scalarField& psi_s,
40  const scalarField& source,
41  const direction cmpt
42 ) const
43 {
45  solveScalarField& psi = tpsi.ref();
46 
48 
49  // Setup class containing solver performance data
50  solverPerformance solverPerf(typeName, fieldName_);
51 
52  // Calculate A.psi used to calculate the initial residual
53  solveScalarField Apsi(psi.size());
55 
56  // Create the storage for the finestCorrection which may be used as a
57  // temporary in normFactor
58  solveScalarField finestCorrection(psi.size());
59 
60  // Calculate normalisation factor
61  solveScalar normFactor =
62  this->normFactor(psi, tsource(), Apsi, finestCorrection);
63 
64  if ((log_ >= 2) || (debug >= 2))
65  {
66  Pout<< " Normalisation factor = " << normFactor << endl;
67  }
68 
69  // Calculate initial finest-grid residual field
70  solveScalarField finestResidual(tsource() - Apsi);
71 
73  (
75  fieldName_,
76  true
77  );
78 
79  // Calculate normalised residual for convergence test
80  solverPerf.initialResidual() = gSumMag
81  (
82  finestResidual,
83  matrix().mesh().comm()
84  )/normFactor;
85  solverPerf.finalResidual() = solverPerf.initialResidual();
86 
87 
88  // Check convergence, solve if not converged
89  if
90  (
91  minIter_ > 0
92  || !solverPerf.checkConvergence(tolerance_, relTol_, log_)
93  )
94  {
95  // Create coarse grid correction fields
96  PtrList<solveScalarField> coarseCorrFields;
97 
98  // Create coarse grid sources
99  PtrList<solveScalarField> coarseSources;
100 
101  // Create the smoothers for all levels
103 
104  // Scratch fields if processor-agglomerated coarse level meshes
105  // are bigger than original. Usually not needed
106  solveScalarField scratch1;
107  solveScalarField scratch2;
108 
109  // Initialise the above data structures
110  initVcycle
111  (
112  coarseCorrFields,
113  coarseSources,
114  smoothers,
115  scratch1,
116  scratch2
117  );
118 
119  do
120  {
121  Vcycle
122  (
123  smoothers,
124  psi,
125  source,
126  Apsi,
127  finestCorrection,
128  finestResidual,
129 
130  (scratch1.size() ? scratch1 : Apsi),
131  (scratch2.size() ? scratch2 : finestCorrection),
132 
133  coarseCorrFields,
134  coarseSources,
135  cmpt
136  );
137 
138  // Calculate finest level residual field
140  finestResidual = tsource();
141  finestResidual -= Apsi;
142 
143  solverPerf.finalResidual() = gSumMag
144  (
145  finestResidual,
146  matrix().mesh().comm()
147  )/normFactor;
148 
149  if ((log_ >= 2) || (debug >= 2))
150  {
151  solverPerf.print(Info.masterStream(matrix().mesh().comm()));
152  }
153  } while
154  (
155  (
156  ++solverPerf.nIterations() < maxIter_
157  && !solverPerf.checkConvergence(tolerance_, relTol_, log_)
158  )
159  || solverPerf.nIterations() < minIter_
160  );
161  }
162 
164  (
166  fieldName_,
167  false
168  );
169 
170  return solverPerf;
171 }
172 
173 
174 void Foam::GAMGSolver::Vcycle
175 (
176  const PtrList<lduMatrix::smoother>& smoothers,
178  const scalarField& source,
179  solveScalarField& Apsi,
180  solveScalarField& finestCorrection,
181  solveScalarField& finestResidual,
182 
183  solveScalarField& scratch1,
184  solveScalarField& scratch2,
185 
186  PtrList<solveScalarField>& coarseCorrFields,
187  PtrList<solveScalarField>& coarseSources,
188  const direction cmpt
189 ) const
190 {
191  //debug = 2;
192 
193  const label coarsestLevel = matrixLevels_.size() - 1;
194 
195  // Restrict finest grid residual for the next level up.
196  agglomeration_.restrictField(coarseSources[0], finestResidual, 0, true);
197 
198  if (nPreSweeps_ && ((log_ >= 2) || (debug >= 2)))
199  {
200  Pout<< "Pre-smoothing scaling factors: ";
201  }
202 
203 
204  // Residual restriction (going to coarser levels)
205  for (label leveli = 0; leveli < coarsestLevel; leveli++)
206  {
207  if (coarseSources.set(leveli + 1))
208  {
209  // If the optional pre-smoothing sweeps are selected
210  // smooth the coarse-grid field for the restricted source
211  if (nPreSweeps_)
212  {
213  coarseCorrFields[leveli] = 0.0;
214 
215  smoothers[leveli + 1].scalarSmooth
216  (
217  coarseCorrFields[leveli],
218  coarseSources[leveli], //coarseSource,
219  cmpt,
220  min
221  (
222  nPreSweeps_ + preSweepsLevelMultiplier_*leveli,
223  maxPreSweeps_
224  )
225  );
226 
227  // Scale coarse-grid correction field
228  // but not on the coarsest level because it evaluates to 1
229  if (scaleCorrection_ && leveli < coarsestLevel - 1)
230  {
232  (
233  scratch1,
234  coarseCorrFields[leveli].size()
235  );
236 
237  scale
238  (
239  coarseCorrFields[leveli],
240  const_cast<solveScalarField&>
241  (
242  ACf.operator const solveScalarField&()
243  ),
244  matrixLevels_[leveli],
245  interfaceLevelsBouCoeffs_[leveli],
246  interfaceLevels_[leveli],
247  coarseSources[leveli],
248  cmpt
249  );
250  }
251 
252  // Correct the residual with the new solution
253  // residual can be used by fusing Amul with b-Amul
254  matrixLevels_[leveli].residual
255  (
256  coarseSources[leveli],
257  coarseCorrFields[leveli],
258  ConstPrecisionAdaptor<scalar, solveScalar>
259  (
260  coarseSources[leveli]
261  )(),
262  interfaceLevelsBouCoeffs_[leveli],
263  interfaceLevels_[leveli],
264  cmpt
265  );
266  }
267 
268  // Residual is equal to source
269  agglomeration_.restrictField
270  (
271  coarseSources[leveli + 1],
272  coarseSources[leveli],
273  leveli + 1,
274  true
275  );
276  }
277  }
278 
279  if (nPreSweeps_ && ((log_ >= 2) || (debug >= 2)))
280  {
281  Pout<< endl;
282  }
283 
284 
285  // Solve Coarsest level with either an iterative or direct solver
286  if (coarseCorrFields.set(coarsestLevel))
287  {
288  solveCoarsestLevel
289  (
290  coarseCorrFields[coarsestLevel],
291  coarseSources[coarsestLevel]
292  );
293  }
294 
295  if ((log_ >= 2) || (debug >= 2))
296  {
297  Pout<< "Post-smoothing scaling factors: ";
298  }
299 
300  // Smoothing and prolongation of the coarse correction fields
301  // (going to finer levels)
302 
303  solveScalarField dummyField(0);
304 
305  // Work storage for prolongation
306  solveScalarField work;
307 
308  for (label leveli = coarsestLevel - 1; leveli >= 0; leveli--)
309  {
310  if (coarseCorrFields.set(leveli))
311  {
312  // Create a field for the pre-smoothed correction field
313  // as a sub-field of the finestCorrection which is not
314  // currently being used
315  solveScalarField::subField preSmoothedCoarseCorrField
316  (
317  scratch2,
318  coarseCorrFields[leveli].size()
319  );
320 
321  // Only store the preSmoothedCoarseCorrField if pre-smoothing is
322  // used
323  if (nPreSweeps_)
324  {
325  preSmoothedCoarseCorrField = coarseCorrFields[leveli];
326  }
327 
328 
329  // Prolong correction to leveli
330  const auto& cf = agglomeration_.prolongField
331  (
332  coarseCorrFields[leveli], // current level
333  work,
334  (
335  coarseCorrFields.set(leveli + 1)
336  ? coarseCorrFields[leveli + 1]
337  : dummyField // dummy value
338  ),
339  leveli + 1
340  );
341 
342 
343  // Create A.psi for this coarse level as a sub-field of Apsi
345  (
346  scratch1,
347  coarseCorrFields[leveli].size()
348  );
349  solveScalarField& ACfRef =
350  const_cast<solveScalarField&>
351  (
352  ACf.operator const solveScalarField&()
353  );
354 
355  if (interpolateCorrection_)
356  {
357  // Normal operation : have both coarse level and fine
358  // level. No processor agglomeration
360  (
361  coarseCorrFields[leveli],
362  ACfRef,
363  matrixLevels_[leveli],
364  interfaceLevelsBouCoeffs_[leveli],
365  interfaceLevels_[leveli],
366  agglomeration_.restrictAddressing(leveli + 1),
367  cf,
368  cmpt
369  );
370  }
371 
372  // Scale coarse-grid correction field
373  // but not on the coarsest level because it evaluates to 1
374  if
375  (
376  scaleCorrection_
377  && (interpolateCorrection_ || leveli < coarsestLevel - 1)
378  )
379  {
380  scale
381  (
382  coarseCorrFields[leveli],
383  ACfRef,
384  matrixLevels_[leveli],
385  interfaceLevelsBouCoeffs_[leveli],
386  interfaceLevels_[leveli],
387  coarseSources[leveli],
388  cmpt
389  );
390  }
391 
392  // Only add the preSmoothedCoarseCorrField if pre-smoothing is
393  // used
394  if (nPreSweeps_)
395  {
396  coarseCorrFields[leveli] += preSmoothedCoarseCorrField;
397  }
398 
399  smoothers[leveli + 1].scalarSmooth
400  (
401  coarseCorrFields[leveli],
402  coarseSources[leveli], //coarseSource,
403  cmpt,
404  min
405  (
406  nPostSweeps_ + postSweepsLevelMultiplier_*leveli,
407  maxPostSweeps_
408  )
409  );
410  }
411  }
412 
413  // Prolong the finest level correction
414  agglomeration_.prolongField
415  (
416  finestCorrection,
417  coarseCorrFields[0],
418  0,
419  true
420  );
421 
422  if (interpolateCorrection_)
423  {
425  (
426  finestCorrection,
427  Apsi,
428  matrix_,
429  interfaceBouCoeffs_,
430  interfaces_,
431  agglomeration_.restrictAddressing(0),
432  coarseCorrFields[0],
433  cmpt
434  );
435  }
436 
437  if (scaleCorrection_)
438  {
439  // Scale the finest level correction
440  scale
441  (
442  finestCorrection,
443  Apsi,
444  matrix_,
445  interfaceBouCoeffs_,
446  interfaces_,
447  finestResidual,
448  cmpt
449  );
450  }
451 
452  forAll(psi, i)
453  {
454  psi[i] += finestCorrection[i];
455  }
456 
457  smoothers[0].smooth
458  (
459  psi,
460  source,
461  cmpt,
462  nFinestSweeps_
463  );
464 }
465 
466 
467 void Foam::GAMGSolver::initVcycle
468 (
469  PtrList<solveScalarField>& coarseCorrFields,
470  PtrList<solveScalarField>& coarseSources,
471  PtrList<lduMatrix::smoother>& smoothers,
472  solveScalarField& scratch1,
473  solveScalarField& scratch2
474 ) const
475 {
476  label maxSize = matrix_.diag().size();
477 
478  coarseCorrFields.setSize(matrixLevels_.size());
479  coarseSources.setSize(matrixLevels_.size());
480  smoothers.setSize(matrixLevels_.size() + 1);
481 
482  // Create the smoother for the finest level
483  smoothers.set
484  (
485  0,
487  (
488  fieldName_,
489  matrix_,
490  interfaceBouCoeffs_,
491  interfaceIntCoeffs_,
492  interfaces_,
493  controlDict_
494  )
495  );
496 
497  forAll(matrixLevels_, leveli)
498  {
499  if (agglomeration_.nCells(leveli) >= 0)
500  {
501  label nCoarseCells = agglomeration_.nCells(leveli);
502 
503  coarseSources.set(leveli, new solveScalarField(nCoarseCells));
504  }
505 
506  if (matrixLevels_.set(leveli))
507  {
508  const lduMatrix& mat = matrixLevels_[leveli];
509 
510  label nCoarseCells = mat.diag().size();
511 
512  maxSize = max(maxSize, nCoarseCells);
513 
514  coarseCorrFields.set(leveli, new solveScalarField(nCoarseCells));
515 
516  smoothers.set
517  (
518  leveli + 1,
520  (
521  fieldName_,
522  matrixLevels_[leveli],
523  interfaceLevelsBouCoeffs_[leveli],
524  interfaceLevelsIntCoeffs_[leveli],
525  interfaceLevels_[leveli],
526  controlDict_
527  )
528  );
529  }
530  }
531 
532  if (maxSize > matrix_.diag().size())
533  {
534  // Allocate some scratch storage
535  scratch1.setSize(maxSize);
536  scratch2.setSize(maxSize);
537  }
538 }
539 
540 
541 Foam::dictionary Foam::GAMGSolver::PCGsolverDict
542 (
543  const scalar tol,
544  const scalar relTol
545 ) const
546 {
547  dictionary dict(IStringStream("solver PCG; preconditioner DIC;")());
548  dict.add("tolerance", tol);
549  dict.add("relTol", relTol);
550 
551  return dict;
552 }
553 
554 
555 Foam::dictionary Foam::GAMGSolver::PBiCGStabSolverDict
556 (
557  const scalar tol,
558  const scalar relTol
559 ) const
560 {
561  dictionary dict(IStringStream("solver PBiCGStab; preconditioner DILU;")());
562  dict.add("tolerance", tol);
563  dict.add("relTol", relTol);
564 
565  return dict;
566 }
567 
568 
569 void Foam::GAMGSolver::solveCoarsestLevel
570 (
571  solveScalarField& coarsestCorrField,
572  const solveScalarField& coarsestSource
573 ) const
574 {
575  const label coarsestLevel = matrixLevels_.size() - 1;
576 
577  const label coarseComm = matrixLevels_[coarsestLevel].mesh().comm();
578 
579  if (directSolveCoarsest_)
580  {
581  PrecisionAdaptor<scalar, solveScalar> tcorrField(coarsestCorrField);
582 
583  coarsestLUMatrixPtr_->solve
584  (
585  tcorrField.ref(),
586  ConstPrecisionAdaptor<scalar, solveScalar>(coarsestSource)()
587  );
588  }
589  //else if
590  //(
591  // agglomeration_.processorAgglomerate()
592  // && procMatrixLevels_.set(coarsestLevel)
593  //)
594  //{
595  // //const labelList& agglomProcIDs = agglomeration_.agglomProcIDs
596  // //(
597  // // coarsestLevel
598  // //);
599  // //
600  // //scalarField allSource;
601  // //
602  // //globalIndex cellOffsets;
603  // //if (Pstream::myProcNo(coarseComm) == agglomProcIDs[0])
604  // //{
605  // // cellOffsets.offsets() =
606  // // agglomeration_.cellOffsets(coarsestLevel);
607  // //}
608  // //
609  // //cellOffsets.gather
610  // //(
611  // // coarseComm,
612  // // agglomProcIDs,
613  // // coarsestSource,
614  // // allSource
615  // //);
616  // //
617  // //scalarField allCorrField;
618  // //solverPerformance coarseSolverPerf;
619  //
620  // label solveComm = agglomeration_.procCommunicator(coarsestLevel);
621  //
622  // coarsestCorrField = 0;
623  // solverPerformance coarseSolverPerf;
624  //
625  // if (Pstream::myProcNo(solveComm) != -1)
626  // {
627  // const lduMatrix& allMatrix = procMatrixLevels_[coarsestLevel];
628  //
629  // {
630  // Pout<< "** Master:Solving on comm:" << solveComm
631  // << " with procs:" << UPstream::procID(solveComm) << endl;
632  //
633  // if (allMatrix.asymmetric())
634  // {
635  // coarseSolverPerf = PBiCGStab
636  // (
637  // "coarsestLevelCorr",
638  // allMatrix,
639  // procInterfaceLevelsBouCoeffs_[coarsestLevel],
640  // procInterfaceLevelsIntCoeffs_[coarsestLevel],
641  // procInterfaceLevels_[coarsestLevel],
642  // PBiCGStabSolverDict(tolerance_, relTol_)
643  // ).solve
644  // (
645  // coarsestCorrField,
646  // coarsestSource
647  // );
648  // }
649  // else
650  // {
651  // coarseSolverPerf = PCG
652  // (
653  // "coarsestLevelCorr",
654  // allMatrix,
655  // procInterfaceLevelsBouCoeffs_[coarsestLevel],
656  // procInterfaceLevelsIntCoeffs_[coarsestLevel],
657  // procInterfaceLevels_[coarsestLevel],
658  // PCGsolverDict(tolerance_, relTol_)
659  // ).solve
660  // (
661  // coarsestCorrField,
662  // coarsestSource
663  // );
664  // }
665  // }
666  // }
667  //
668  // Pout<< "done master solve." << endl;
669  //
670  // //// Scatter to all processors
671  // //coarsestCorrField.setSize(coarsestSource.size());
672  // //cellOffsets.scatter
673  // //(
674  // // coarseComm,
675  // // agglomProcIDs,
676  // // allCorrField,
677  // // coarsestCorrField
678  // //);
679  //
680  // if (debug >= 2)
681  // {
682  // coarseSolverPerf.print(Info.masterStream(coarseComm));
683  // }
684  //
685  // Pout<< "procAgglom: coarsestSource :" << coarsestSource << endl;
686  // Pout<< "procAgglom: coarsestCorrField:" << coarsestCorrField << endl;
687  //}
688  else
689  {
690  coarsestCorrField = 0;
691  const solverPerformance coarseSolverPerf
692  (
693  coarsestSolverPtr_->scalarSolve
694  (
695  coarsestCorrField,
696  coarsestSource
697  )
698  );
699 
700  if ((log_ >= 2) || debug)
701  {
702  coarseSolverPerf.print(Info.masterStream(coarseComm));
703  }
704  }
705 }
706 
707 
708 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:46
Field< solveScalar > solveScalarField
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
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
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
const FieldField< Field, scalar > & interfaceBouCoeffs_
Definition: lduMatrix.H:137
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:625
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve.
A non-const Field/List wrapper with possible data conversion.
bool checkConvergence(const Type &tolerance, const Type &relTolerance, const int logLevel=0)
Check, store and return convergence.
int log_
Verbosity level for solver output statements.
Definition: lduMatrix.H:149
label minIter_
Minimum number of iterations in the solver.
Definition: lduMatrix.H:154
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:169
const labelType & nIterations() const noexcept
Return number of iterations.
void Amul(solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
Matrix multiplication with updated interfaces.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
label maxIter_
Maximum number of iterations in the solver.
Definition: lduMatrix.H:159
A const Field/List wrapper with possible data conversion.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
const Type & finalResidual() const noexcept
Return final residual.
const lduMatrix & matrix_
Definition: lduMatrix.H:136
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
lduInterfaceFieldPtrsList interfaces_
Definition: lduMatrix.H:139
solveScalarField::cmptType normFactor(const solveScalarField &psi, const solveScalarField &source, const solveScalarField &Apsi, solveScalarField &tmpField, const lduMatrix::normTypes normType) const
Return the matrix norm using the specified norm method.
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
void print(Ostream &os) const
Print summary of solver performance to the given stream.
OSstream & masterStream(const label communicator)
Return OSstream for output operations on the master process only, Snull on other processes.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
static autoPtr< smoother > New(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Return a new smoother.
messageStream Info
Information stream (stdout output on master, null elsewhere)
SubField< solveScalar > subField
Declare type of subField.
Definition: Field.H:128
const volScalarField & psi
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
scalar relTol_
Convergence tolerance relative to the initial.
Definition: lduMatrix.H:174
const Type & initialResidual() const noexcept
Return initial residual.