GAMGSolver.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) 2021-2022 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 "GAMGSolver.H"
30 #include "GAMGInterface.H"
31 #include "PCG.H"
32 #include "PBiCGStab.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
39 
40  lduMatrix::solver::addsymMatrixConstructorToTable<GAMGSolver>
42 
43  lduMatrix::solver::addasymMatrixConstructorToTable<GAMGSolver>
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const word& fieldName,
53  const lduMatrix& matrix,
54  const FieldField<Field, scalar>& interfaceBouCoeffs,
55  const FieldField<Field, scalar>& interfaceIntCoeffs,
56  const lduInterfaceFieldPtrsList& interfaces,
57  const dictionary& solverControls
58 )
59 :
61  (
62  fieldName,
63  matrix,
64  interfaceBouCoeffs,
65  interfaceIntCoeffs,
66  interfaces,
67  solverControls
68  ),
69 
70  // Default values for all controls
71  // which may be overridden by those in controlDict
72  nPreSweeps_(0),
73  preSweepsLevelMultiplier_(1),
74  maxPreSweeps_(4),
75  nPostSweeps_(2),
76  postSweepsLevelMultiplier_(1),
77  maxPostSweeps_(4),
78  nFinestSweeps_(2),
79 
80  cacheAgglomeration_(true),
81  interpolateCorrection_(false),
82  scaleCorrection_(matrix.symmetric()),
83  directSolveCoarsest_(false),
84 
85  agglomeration_(GAMGAgglomeration::New(matrix_, controlDict_)),
86 
87  matrixLevels_(agglomeration_.size()),
88  primitiveInterfaceLevels_(agglomeration_.size()),
89  interfaceLevels_(agglomeration_.size()),
90  interfaceLevelsBouCoeffs_(agglomeration_.size()),
91  interfaceLevelsIntCoeffs_(agglomeration_.size())
92 {
93  readControls();
94 
95  if (agglomeration_.processorAgglomerate())
96  {
97  forAll(agglomeration_, fineLevelIndex)
98  {
99  if (agglomeration_.hasMeshLevel(fineLevelIndex))
100  {
101  if
102  (
103  (fineLevelIndex+1) < agglomeration_.size()
104  && agglomeration_.hasProcMesh(fineLevelIndex+1)
105  )
106  {
107  // Construct matrix without referencing the coarse mesh so
108  // construct a dummy mesh instead. This will get overwritten
109  // by the call to procAgglomerateMatrix so is only to get
110  // it through agglomerateMatrix
111 
112 
113  const lduInterfacePtrsList& fineMeshInterfaces =
114  agglomeration_.interfaceLevel(fineLevelIndex);
115 
116  PtrList<GAMGInterface> dummyPrimMeshInterfaces
117  (
118  fineMeshInterfaces.size()
119  );
120  lduInterfacePtrsList dummyMeshInterfaces
121  (
122  dummyPrimMeshInterfaces.size()
123  );
124  forAll(fineMeshInterfaces, intI)
125  {
126  if (fineMeshInterfaces.set(intI))
127  {
129  refCast<const GAMGInterface>
130  (
131  fineMeshInterfaces[intI]
132  ).write(os);
134 
135  dummyPrimMeshInterfaces.set
136  (
137  intI,
139  (
140  fineMeshInterfaces[intI].type(),
141  intI,
142  dummyMeshInterfaces,
143  is
144  )
145  );
146  }
147  }
148 
149  forAll(dummyPrimMeshInterfaces, intI)
150  {
151  if (dummyPrimMeshInterfaces.set(intI))
152  {
153  dummyMeshInterfaces.set
154  (
155  intI,
156  &dummyPrimMeshInterfaces[intI]
157  );
158  }
159  }
160 
161  // So:
162  // - pass in incorrect mesh (= fine mesh instead of coarse)
163  // - pass in dummy interfaces
164  agglomerateMatrix
165  (
166  fineLevelIndex,
167  agglomeration_.meshLevel(fineLevelIndex),
168  dummyMeshInterfaces
169  );
170 
171 
172  const labelList& procAgglomMap =
173  agglomeration_.procAgglomMap(fineLevelIndex+1);
174  const List<label>& procIDs =
175  agglomeration_.agglomProcIDs(fineLevelIndex+1);
176 
177  procAgglomerateMatrix
178  (
179  procAgglomMap,
180  procIDs,
181  fineLevelIndex
182  );
183  }
184  else
185  {
186  agglomerateMatrix
187  (
188  fineLevelIndex,
189  agglomeration_.meshLevel(fineLevelIndex + 1),
190  agglomeration_.interfaceLevel(fineLevelIndex + 1)
191  );
192  }
193  }
194  else
195  {
196  // No mesh. Not involved in calculation anymore
197  }
198  }
199  }
200  else
201  {
202  forAll(agglomeration_, fineLevelIndex)
203  {
204  // Agglomerate on to coarse level mesh
205  agglomerateMatrix
206  (
207  fineLevelIndex,
208  agglomeration_.meshLevel(fineLevelIndex + 1),
209  agglomeration_.interfaceLevel(fineLevelIndex + 1)
210  );
211  }
212  }
213 
214  if ((log_ >= 2) || (debug & 2))
215  {
216  for
217  (
218  label fineLevelIndex = 0;
219  fineLevelIndex <= matrixLevels_.size();
220  fineLevelIndex++
221  )
222  {
223  if (fineLevelIndex == 0 || matrixLevels_.set(fineLevelIndex-1))
224  {
225  const lduMatrix& matrix = matrixLevel(fineLevelIndex);
227  interfaceLevel(fineLevelIndex);
228 
229  Pout<< "level:" << fineLevelIndex << nl
230  << " nCells:" << matrix.diag().size() << nl
231  << " nFaces:" << matrix.lower().size() << nl
232  << " nInterfaces:" << interfaces.size()
233  << endl;
234 
235  forAll(interfaces, i)
236  {
237  if (interfaces.set(i))
238  {
239  Pout<< " " << i
240  << "\ttype:" << interfaces[i].type()
241  << "\tsize:"
242  << interfaces[i].interface().faceCells().size()
243  << endl;
244  }
245  }
246  }
247  else
248  {
249  Pout<< "level:" << fineLevelIndex << " : no matrix" << endl;
250  }
251  }
252  Pout<< endl;
253  }
254 
255 
256  if (matrixLevels_.size())
257  {
258  const label coarsestLevel = matrixLevels_.size() - 1;
259 
260  if (matrixLevels_.set(coarsestLevel))
261  {
262  if (directSolveCoarsest_)
263  {
264  coarsestLUMatrixPtr_.reset
265  (
266  new LUscalarMatrix
267  (
268  matrixLevels_[coarsestLevel],
269  interfaceLevelsBouCoeffs_[coarsestLevel],
270  interfaceLevels_[coarsestLevel]
271  )
272  );
273  }
274  else
275  {
276  entry* coarseEntry = controlDict_.findEntry
277  (
278  "coarsestLevelCorr",
280  );
281  if (coarseEntry && coarseEntry->isDict())
282  {
283  coarsestSolverPtr_ = lduMatrix::solver::New
284  (
285  "coarsestLevelCorr",
286  matrixLevels_[coarsestLevel],
287  interfaceLevelsBouCoeffs_[coarsestLevel],
288  interfaceLevelsIntCoeffs_[coarsestLevel],
289  interfaceLevels_[coarsestLevel],
290  coarseEntry->dict()
291  );
292  }
293  else if (matrixLevels_[coarsestLevel].asymmetric())
294  {
295  coarsestSolverPtr_.reset
296  (
297  new PBiCGStab
298  (
299  "coarsestLevelCorr",
300  matrixLevels_[coarsestLevel],
301  interfaceLevelsBouCoeffs_[coarsestLevel],
302  interfaceLevelsIntCoeffs_[coarsestLevel],
303  interfaceLevels_[coarsestLevel],
304  PBiCGStabSolverDict(tolerance_, relTol_)
305  )
306  );
307  }
308  else
309  {
310  coarsestSolverPtr_.reset
311  (
312  new PCG
313  (
314  "coarsestLevelCorr",
315  matrixLevels_[coarsestLevel],
316  interfaceLevelsBouCoeffs_[coarsestLevel],
317  interfaceLevelsIntCoeffs_[coarsestLevel],
318  interfaceLevels_[coarsestLevel],
319  PCGsolverDict(tolerance_, relTol_)
320  )
321  );
322  }
323  }
324  }
325  }
326  else
327  {
329  << "No coarse levels created, either matrix too small for GAMG"
330  " or nCellsInCoarsestLevel too large.\n"
331  " Either choose another solver of reduce "
332  "nCellsInCoarsestLevel."
334  }
335 }
336 
337 
338 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
339 
341 {
342  if (!cacheAgglomeration_)
343  {
344  delete &agglomeration_;
345  }
346 }
347 
348 
349 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
350 
351 void Foam::GAMGSolver::readControls()
352 {
354 
355  controlDict_.readIfPresent("cacheAgglomeration", cacheAgglomeration_);
356  controlDict_.readIfPresent("nPreSweeps", nPreSweeps_);
357  controlDict_.readIfPresent
358  (
359  "preSweepsLevelMultiplier",
360  preSweepsLevelMultiplier_
361  );
362  controlDict_.readIfPresent("maxPreSweeps", maxPreSweeps_);
363  controlDict_.readIfPresent("nPostSweeps", nPostSweeps_);
364  controlDict_.readIfPresent
365  (
366  "postSweepsLevelMultiplier",
367  postSweepsLevelMultiplier_
368  );
369  controlDict_.readIfPresent("maxPostSweeps", maxPostSweeps_);
370  controlDict_.readIfPresent("nFinestSweeps", nFinestSweeps_);
371  controlDict_.readIfPresent("interpolateCorrection", interpolateCorrection_);
372  controlDict_.readIfPresent("scaleCorrection", scaleCorrection_);
373  controlDict_.readIfPresent("directSolveCoarsest", directSolveCoarsest_);
374 
375  if ((log_ >= 2) || debug)
376  {
377  Info<< "GAMGSolver settings :"
378  << " cacheAgglomeration:" << cacheAgglomeration_
379  << " nPreSweeps:" << nPreSweeps_
380  << " preSweepsLevelMultiplier:" << preSweepsLevelMultiplier_
381  << " maxPreSweeps:" << maxPreSweeps_
382  << " nPostSweeps:" << nPostSweeps_
383  << " postSweepsLevelMultiplier:" << postSweepsLevelMultiplier_
384  << " maxPostSweeps:" << maxPostSweeps_
385  << " nFinestSweeps:" << nFinestSweeps_
386  << " interpolateCorrection:" << interpolateCorrection_
387  << " scaleCorrection:" << scaleCorrection_
388  << " directSolveCoarsest:" << directSolveCoarsest_
389  << endl;
390  }
391 }
392 
393 
394 const Foam::lduMatrix& Foam::GAMGSolver::matrixLevel(const label i) const
395 {
396  return i ? matrixLevels_[i-1] : matrix_;
397 }
398 
399 
400 const Foam::lduInterfaceFieldPtrsList& Foam::GAMGSolver::interfaceLevel
401 (
402  const label i
403 ) const
404 {
405  return i ? interfaceLevels_[i-1] : interfaces_;
406 }
407 
408 
410 Foam::GAMGSolver::interfaceBouCoeffsLevel
411 (
412  const label i
413 ) const
414 {
415  return i ? interfaceLevelsBouCoeffs_[i-1] : interfaceBouCoeffs_;
416 }
417 
418 
420 Foam::GAMGSolver::interfaceIntCoeffsLevel
421 (
422  const label i
423 ) const
424 {
425  return i ? interfaceLevelsIntCoeffs_[i-1] : interfaceIntCoeffs_;
426 }
427 
428 
429 // ************************************************************************* //
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
bool hasProcMesh(const label fineLeveli) const
Check that level has combined mesh.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
static autoPtr< GAMGInterface > New(const label index, const lduInterfacePtrsList &coarseInterfaces, const lduInterface &fineInterface, const labelField &localRestrictAddressing, const labelField &neighbourRestrictAddressing, const label fineLevelIndex, const label coarseComm)
Return a pointer to a new interface created on freestore given.
GAMGSolver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Construct from lduMatrix and solver controls.
Definition: GAMGSolver.C:44
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Base solver class.
Definition: solver.H:45
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const lduMesh & meshLevel(const label leveli) const
Return LDU mesh of given level.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
int log_
Verbosity level for solver output statements.
Definition: lduMatrix.H:149
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:169
bool processorAgglomerate() const
Whether to agglomerate across processors.
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const labelList & agglomProcIDs(const label fineLeveli) const
Set of processors to agglomerate. Element 0 is the master processor. (local, same only on those proce...
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:321
virtual const word & type() const =0
Runtime type information.
lduMatrix::solver::addsymMatrixConstructorToTable< GAMGSolver > addGAMGSolverMatrixConstructorToTable_
Definition: GAMGSolver.C:34
static autoPtr< solver > New(const word &solverName, 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 solver of given type.
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
dictionary controlDict_
Dictionary of solution controls.
Definition: lduMatrix.H:144
Geometric agglomerated algebraic multigrid solver.
Definition: GAMGSolver.H:72
virtual ~GAMGSolver()
Destructor.
Definition: GAMGSolver.C:333
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
scalarField & lower()
Definition: lduMatrix.C:179
const labelList & procAgglomMap(const label fineLeveli) const
Mapping from processor to agglomerated processor (global, all processors have the same information)...
Input from string buffer, using a ISstream. Always UNCOMPRESSED.
Definition: StringStream.H:120
lduMatrix::solver::addasymMatrixConstructorToTable< GAMGSolver > addGAMGAsymSolverMatrixConstructorToTable_
Definition: GAMGSolver.C:37
const lduInterfacePtrsList & interfaceLevel(const label leveli) const
Return LDU interface addressing of given level.
const entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:84
messageStream Info
Information stream (stdout output on master, null elsewhere)
scalarField & diag()
Definition: lduMatrix.C:197
Geometric agglomerated algebraic multigrid agglomeration class.
Output to string buffer, using a OSstream. Always UNCOMPRESSED.
Definition: StringStream.H:256
bool hasMeshLevel(const label leveli) const
Do we have mesh for given level?
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
virtual void readControls()
Read the control parameters from controlDict_.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: UPtrList.H:366
scalar relTol_
Convergence tolerance relative to the initial.
Definition: lduMatrix.H:174