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 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
36namespace 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 cacheAgglomeration_(true),
73 nPreSweeps_(0),
74 preSweepsLevelMultiplier_(1),
75 maxPreSweeps_(4),
76 nPostSweeps_(2),
77 postSweepsLevelMultiplier_(1),
78 maxPostSweeps_(4),
79 nFinestSweeps_(2),
80 interpolateCorrection_(false),
81 scaleCorrection_(matrix.symmetric()),
82 directSolveCoarsest_(false),
83 agglomeration_(GAMGAgglomeration::New(matrix_, controlDict_)),
84
85 matrixLevels_(agglomeration_.size()),
86 primitiveInterfaceLevels_(agglomeration_.size()),
87 interfaceLevels_(agglomeration_.size()),
88 interfaceLevelsBouCoeffs_(agglomeration_.size()),
89 interfaceLevelsIntCoeffs_(agglomeration_.size())
90{
91 readControls();
92
93 if (agglomeration_.processorAgglomerate())
94 {
95 forAll(agglomeration_, fineLevelIndex)
96 {
97 if (agglomeration_.hasMeshLevel(fineLevelIndex))
98 {
99 if
100 (
101 (fineLevelIndex+1) < agglomeration_.size()
102 && agglomeration_.hasProcMesh(fineLevelIndex+1)
103 )
104 {
105 // Construct matrix without referencing the coarse mesh so
106 // construct a dummy mesh instead. This will get overwritten
107 // by the call to procAgglomerateMatrix so is only to get
108 // it through agglomerateMatrix
109
110
111 const lduInterfacePtrsList& fineMeshInterfaces =
112 agglomeration_.interfaceLevel(fineLevelIndex);
113
114 PtrList<GAMGInterface> dummyPrimMeshInterfaces
115 (
116 fineMeshInterfaces.size()
117 );
118 lduInterfacePtrsList dummyMeshInterfaces
119 (
120 dummyPrimMeshInterfaces.size()
121 );
122 forAll(fineMeshInterfaces, intI)
123 {
124 if (fineMeshInterfaces.set(intI))
125 {
127 refCast<const GAMGInterface>
128 (
129 fineMeshInterfaces[intI]
130 ).write(os);
131 IStringStream is(os.str());
132
133 dummyPrimMeshInterfaces.set
134 (
135 intI,
137 (
138 fineMeshInterfaces[intI].type(),
139 intI,
140 dummyMeshInterfaces,
141 is
142 )
143 );
144 }
145 }
146
147 forAll(dummyPrimMeshInterfaces, intI)
148 {
149 if (dummyPrimMeshInterfaces.set(intI))
150 {
151 dummyMeshInterfaces.set
152 (
153 intI,
154 &dummyPrimMeshInterfaces[intI]
155 );
156 }
157 }
158
159 // So:
160 // - pass in incorrect mesh (= fine mesh instead of coarse)
161 // - pass in dummy interfaces
162 agglomerateMatrix
163 (
164 fineLevelIndex,
165 agglomeration_.meshLevel(fineLevelIndex),
166 dummyMeshInterfaces
167 );
168
169
170 const labelList& procAgglomMap =
171 agglomeration_.procAgglomMap(fineLevelIndex+1);
172 const List<label>& procIDs =
173 agglomeration_.agglomProcIDs(fineLevelIndex+1);
174
175 procAgglomerateMatrix
176 (
177 procAgglomMap,
178 procIDs,
179 fineLevelIndex
180 );
181 }
182 else
183 {
184 agglomerateMatrix
185 (
186 fineLevelIndex,
187 agglomeration_.meshLevel(fineLevelIndex + 1),
188 agglomeration_.interfaceLevel(fineLevelIndex + 1)
189 );
190 }
191 }
192 else
193 {
194 // No mesh. Not involved in calculation anymore
195 }
196 }
197 }
198 else
199 {
200 forAll(agglomeration_, fineLevelIndex)
201 {
202 // Agglomerate on to coarse level mesh
203 agglomerateMatrix
204 (
205 fineLevelIndex,
206 agglomeration_.meshLevel(fineLevelIndex + 1),
207 agglomeration_.interfaceLevel(fineLevelIndex + 1)
208 );
209 }
210 }
211
212 if ((log_ >= 2) || (debug & 2))
213 {
214 for
215 (
216 label fineLevelIndex = 0;
217 fineLevelIndex <= matrixLevels_.size();
218 fineLevelIndex++
219 )
220 {
221 if (fineLevelIndex == 0 || matrixLevels_.set(fineLevelIndex-1))
222 {
223 const lduMatrix& matrix = matrixLevel(fineLevelIndex);
225 interfaceLevel(fineLevelIndex);
226
227 Pout<< "level:" << fineLevelIndex << nl
228 << " nCells:" << matrix.diag().size() << nl
229 << " nFaces:" << matrix.lower().size() << nl
230 << " nInterfaces:" << interfaces.size()
231 << endl;
232
234 {
235 if (interfaces.set(i))
236 {
237 Pout<< " " << i
238 << "\ttype:" << interfaces[i].type()
239 << endl;
240 }
241 }
242 }
243 else
244 {
245 Pout<< "level:" << fineLevelIndex << " : no matrix" << endl;
246 }
247 }
248 Pout<< endl;
249 }
250
251
252 if (matrixLevels_.size())
253 {
254 const label coarsestLevel = matrixLevels_.size() - 1;
255
256 if (matrixLevels_.set(coarsestLevel))
257 {
258 if (directSolveCoarsest_)
259 {
260 coarsestLUMatrixPtr_.reset
261 (
263 (
264 matrixLevels_[coarsestLevel],
265 interfaceLevelsBouCoeffs_[coarsestLevel],
266 interfaceLevels_[coarsestLevel]
267 )
268 );
269 }
270 else
271 {
272 entry* coarseEntry = controlDict_.findEntry
273 (
274 "coarsestLevelCorr",
276 );
277 if (coarseEntry && coarseEntry->isDict())
278 {
279 coarsestSolverPtr_ = lduMatrix::solver::New
280 (
281 "coarsestLevelCorr",
282 matrixLevels_[coarsestLevel],
283 interfaceLevelsBouCoeffs_[coarsestLevel],
284 interfaceLevelsIntCoeffs_[coarsestLevel],
285 interfaceLevels_[coarsestLevel],
286 coarseEntry->dict()
287 );
288 }
289 else if (matrixLevels_[coarsestLevel].asymmetric())
290 {
291 coarsestSolverPtr_.reset
292 (
293 new PBiCGStab
294 (
295 "coarsestLevelCorr",
296 matrixLevels_[coarsestLevel],
297 interfaceLevelsBouCoeffs_[coarsestLevel],
298 interfaceLevelsIntCoeffs_[coarsestLevel],
299 interfaceLevels_[coarsestLevel],
300 PBiCGStabSolverDict(tolerance_, relTol_)
301 )
302 );
303 }
304 else
305 {
306 coarsestSolverPtr_.reset
307 (
308 new PCG
309 (
310 "coarsestLevelCorr",
311 matrixLevels_[coarsestLevel],
312 interfaceLevelsBouCoeffs_[coarsestLevel],
313 interfaceLevelsIntCoeffs_[coarsestLevel],
314 interfaceLevels_[coarsestLevel],
315 PCGsolverDict(tolerance_, relTol_)
316 )
317 );
318 }
319 }
320 }
321 }
322 else
323 {
325 << "No coarse levels created, either matrix too small for GAMG"
326 " or nCellsInCoarsestLevel too large.\n"
327 " Either choose another solver of reduce "
328 "nCellsInCoarsestLevel."
329 << exit(FatalError);
330 }
331}
332
333
334// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
335
337{
338 if (!cacheAgglomeration_)
339 {
340 delete &agglomeration_;
341 }
342}
343
344
345// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
346
347void Foam::GAMGSolver::readControls()
348{
350
351 controlDict_.readIfPresent("cacheAgglomeration", cacheAgglomeration_);
352 controlDict_.readIfPresent("nPreSweeps", nPreSweeps_);
353 controlDict_.readIfPresent
354 (
355 "preSweepsLevelMultiplier",
356 preSweepsLevelMultiplier_
357 );
358 controlDict_.readIfPresent("maxPreSweeps", maxPreSweeps_);
359 controlDict_.readIfPresent("nPostSweeps", nPostSweeps_);
360 controlDict_.readIfPresent
361 (
362 "postSweepsLevelMultiplier",
363 postSweepsLevelMultiplier_
364 );
365 controlDict_.readIfPresent("maxPostSweeps", maxPostSweeps_);
366 controlDict_.readIfPresent("nFinestSweeps", nFinestSweeps_);
367 controlDict_.readIfPresent("interpolateCorrection", interpolateCorrection_);
368 controlDict_.readIfPresent("scaleCorrection", scaleCorrection_);
369 controlDict_.readIfPresent("directSolveCoarsest", directSolveCoarsest_);
370
371 if ((log_ >= 2) || debug)
372 {
373 Info<< "GAMGSolver settings :"
374 << " cacheAgglomeration:" << cacheAgglomeration_
375 << " nPreSweeps:" << nPreSweeps_
376 << " preSweepsLevelMultiplier:" << preSweepsLevelMultiplier_
377 << " maxPreSweeps:" << maxPreSweeps_
378 << " nPostSweeps:" << nPostSweeps_
379 << " postSweepsLevelMultiplier:" << postSweepsLevelMultiplier_
380 << " maxPostSweeps:" << maxPostSweeps_
381 << " nFinestSweeps:" << nFinestSweeps_
382 << " interpolateCorrection:" << interpolateCorrection_
383 << " scaleCorrection:" << scaleCorrection_
384 << " directSolveCoarsest:" << directSolveCoarsest_
385 << endl;
386 }
387}
388
389
390const Foam::lduMatrix& Foam::GAMGSolver::matrixLevel(const label i) const
391{
392 if (i == 0)
393 {
394 return matrix_;
395 }
396 else
397 {
398 return matrixLevels_[i - 1];
399 }
400}
401
402
403const Foam::lduInterfaceFieldPtrsList& Foam::GAMGSolver::interfaceLevel
404(
405 const label i
406) const
407{
408 if (i == 0)
409 {
410 return interfaces_;
411 }
412 else
413 {
414 return interfaceLevels_[i - 1];
415 }
416}
417
418
420Foam::GAMGSolver::interfaceBouCoeffsLevel
421(
422 const label i
423) const
424{
425 if (i == 0)
426 {
427 return interfaceBouCoeffs_;
428 }
429 else
430 {
431 return interfaceLevelsBouCoeffs_[i - 1];
432 }
433}
434
435
437Foam::GAMGSolver::interfaceIntCoeffsLevel
438(
439 const label i
440) const
441{
442 if (i == 0)
443 {
444 return interfaceIntCoeffs_;
445 }
446 else
447 {
448 return interfaceLevelsIntCoeffs_[i - 1];
449 }
450}
451
452
453// ************************************************************************* //
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
Geometric agglomerated algebraic multigrid agglomeration class.
bool processorAgglomerate() const
Whether to agglomerate across processors.
const labelList & procAgglomMap(const label fineLeveli) const
bool hasProcMesh(const label fineLeveli) const
Check that level has combined mesh.
const labelList & agglomProcIDs(const label fineLeveli) const
const lduInterfacePtrsList & interfaceLevel(const label leveli) const
Return LDU interface addressing of given level.
bool hasMeshLevel(const label leveli) const
Do we have mesh for given level?
const lduMesh & meshLevel(const label leveli) const
Return LDU mesh of given level.
Geometric agglomerated algebraic multigrid solver.
Definition: GAMGSolver.H:79
virtual ~GAMGSolver()
Destructor.
Definition: GAMGSolver.C:336
Input from string buffer, using a ISstream. Always UNCOMPRESSED.
Definition: StringStream.H:112
Class to perform the LU decomposition on a symmetric matrix.
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
Output to string buffer, using a OSstream. Always UNCOMPRESSED.
Definition: StringStream.H:231
Preconditioned bi-conjugate gradient stabilized solver for asymmetric lduMatrices using a run-time se...
Definition: PBiCGStab.H:70
Preconditioned conjugate gradient solver for symmetric lduMatrices using a run-time selectable precon...
Definition: PCG.H:58
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const T * set(const label i) const
Definition: UPtrList.H:248
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find for an entry (non-const access) with the given keyword.
Definition: dictionaryI.H:97
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
virtual bool isDict() const noexcept
Return true if this entry is a dictionary.
Definition: entry.H:233
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary.
@ LITERAL_RECURSIVE
Definition: keyType.H:86
scalar tolerance_
Final convergence tolerance.
Definition: lduMatrix.H:126
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:248
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:233
int log_
Level of verbosity in the solver output statements.
Definition: lduMatrix.H:117
scalar relTol_
Convergence tolerance relative to the initial.
Definition: lduMatrix.H:129
virtual void readControls()
Read the control parameters from the controlDict_.
dictionary controlDict_
Dictionary of controls.
Definition: lduMatrix.H:114
virtual const word & type() const =0
Runtime type information.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:84
scalarField & lower()
Definition: lduMatrix.C:174
scalarField & diag()
Definition: lduMatrix.C:192
Base class for solution control classes.
Definition: solver.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
lduMatrix::solver::addasymMatrixConstructorToTable< GAMGSolver > addGAMGAsymSolverMatrixConstructorToTable_
Definition: GAMGSolver.C:44
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
lduMatrix::solver::addsymMatrixConstructorToTable< GAMGSolver > addGAMGSolverMatrixConstructorToTable_
Definition: GAMGSolver.C:41
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333