Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
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{
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
127
128 forAll(fineMeshInterfaces, intI)
129 {
130 if (fineMeshInterfaces.set(intI))
131 {
132 os.rewind();
133
135 (
136 fineMeshInterfaces[intI]
137 ).write(os);
138
139 is.reset(os.view());
140
141 dummyPrimMeshInterfaces.set
142 (
143 intI,
145 (
146 fineMeshInterfaces[intI].type(),
147 intI,
148 dummyMeshInterfaces,
149 is
150 )
151 );
152 }
153 }
154
155 forAll(dummyPrimMeshInterfaces, intI)
156 {
157 if (dummyPrimMeshInterfaces.set(intI))
158 {
159 dummyMeshInterfaces.set
160 (
161 intI,
162 &dummyPrimMeshInterfaces[intI]
163 );
164 }
165 }
166
167 // So:
168 // - pass in incorrect mesh (= fine mesh instead of coarse)
169 // - pass in dummy interfaces
170 agglomerateMatrix
171 (
172 fineLevelIndex,
173 agglomeration_.meshLevel(fineLevelIndex),
174 dummyMeshInterfaces
175 );
176
177
178 const labelList& procAgglomMap =
179 agglomeration_.procAgglomMap(fineLevelIndex+1);
180 const List<label>& procIDs =
181 agglomeration_.agglomProcIDs(fineLevelIndex+1);
182
183 procAgglomerateMatrix
184 (
185 procAgglomMap,
186 procIDs,
187 fineLevelIndex
188 );
189 }
190 else
191 {
192 agglomerateMatrix
193 (
194 fineLevelIndex,
195 agglomeration_.meshLevel(fineLevelIndex + 1),
196 agglomeration_.interfaceLevel(fineLevelIndex + 1)
197 );
198 }
199 }
200 else
201 {
202 // No mesh. Not involved in calculation anymore
203 }
204 }
205 }
206 else
207 {
208 forAll(agglomeration_, fineLevelIndex)
209 {
210 // Agglomerate on to coarse level mesh
211 agglomerateMatrix
212 (
213 fineLevelIndex,
214 agglomeration_.meshLevel(fineLevelIndex + 1),
215 agglomeration_.interfaceLevel(fineLevelIndex + 1)
216 );
217 }
218 }
219
220 if ((log_ >= 2) || (debug & 2))
221 {
222 for
223 (
224 label fineLevelIndex = 0;
225 fineLevelIndex <= matrixLevels_.size();
226 fineLevelIndex++
227 )
228 {
229 if (fineLevelIndex == 0 || matrixLevels_.set(fineLevelIndex-1))
230 {
231 const lduMatrix& matrix = matrixLevel(fineLevelIndex);
233 interfaceLevel(fineLevelIndex);
234
235 Pout<< "level:" << fineLevelIndex << nl
236 << " nCells:" << matrix.diag().size() << nl
237 << " nFaces:" << matrix.lower().size() << nl
238 << " nInterfaces:" << interfaces.size()
239 << endl;
240
242 {
243 if (interfaces.set(i))
244 {
245 Pout<< " " << i
246 << "\ttype:" << interfaces[i].type()
247 << "\tsize:"
248 << interfaces[i].interface().faceCells().size()
249 << endl;
250 }
251 }
252 }
253 else
254 {
255 Pout<< "level:" << fineLevelIndex << " : no matrix" << endl;
256 }
257 }
258 Pout<< endl;
259 }
260
261
262 if (matrixLevels_.size())
263 {
264 const label coarsestLevel = matrixLevels_.size() - 1;
265
266 if (matrixLevels_.set(coarsestLevel))
267 {
268 if (directSolveCoarsest_)
269 {
270 coarsestLUMatrixPtr_.reset
271 (
272 new LUscalarMatrix
273 (
274 matrixLevels_[coarsestLevel],
275 interfaceLevelsBouCoeffs_[coarsestLevel],
276 interfaceLevels_[coarsestLevel]
277 )
278 );
279 }
280 else
281 {
282 entry* coarseEntry = controlDict_.findEntry
283 (
284 "coarsestLevelCorr",
286 );
287 if (coarseEntry && coarseEntry->isDict())
288 {
289 coarsestSolverPtr_ = lduMatrix::solver::New
290 (
291 "coarsestLevelCorr",
292 matrixLevels_[coarsestLevel],
293 interfaceLevelsBouCoeffs_[coarsestLevel],
294 interfaceLevelsIntCoeffs_[coarsestLevel],
295 interfaceLevels_[coarsestLevel],
296 coarseEntry->dict()
297 );
298 }
299 else if (matrixLevels_[coarsestLevel].asymmetric())
300 {
301 coarsestSolverPtr_.reset
302 (
303 new PBiCGStab
304 (
305 "coarsestLevelCorr",
306 matrixLevels_[coarsestLevel],
307 interfaceLevelsBouCoeffs_[coarsestLevel],
308 interfaceLevelsIntCoeffs_[coarsestLevel],
309 interfaceLevels_[coarsestLevel],
310 PBiCGStabSolverDict(tolerance_, relTol_)
311 )
312 );
313 }
314 else
315 {
316 coarsestSolverPtr_.reset
317 (
318 new PCG
319 (
320 "coarsestLevelCorr",
321 matrixLevels_[coarsestLevel],
322 interfaceLevelsBouCoeffs_[coarsestLevel],
323 interfaceLevelsIntCoeffs_[coarsestLevel],
324 interfaceLevels_[coarsestLevel],
325 PCGsolverDict(tolerance_, relTol_)
326 )
327 );
328 }
329 }
330 }
331 }
332 else
333 {
335 << "No coarse levels created, either matrix too small for GAMG"
336 " or nCellsInCoarsestLevel too large.\n"
337 " Either choose another solver of reduce "
338 "nCellsInCoarsestLevel."
340 }
341}
342
343
344// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
345
347{
348 if (!cacheAgglomeration_)
349 {
350 delete &agglomeration_;
351 }
352}
353
354
355// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
356
358{
360
361 controlDict_.readIfPresent("cacheAgglomeration", cacheAgglomeration_);
362 controlDict_.readIfPresent("nPreSweeps", nPreSweeps_);
363 controlDict_.readIfPresent
364 (
365 "preSweepsLevelMultiplier",
366 preSweepsLevelMultiplier_
367 );
368 controlDict_.readIfPresent("maxPreSweeps", maxPreSweeps_);
369 controlDict_.readIfPresent("nPostSweeps", nPostSweeps_);
370 controlDict_.readIfPresent
371 (
372 "postSweepsLevelMultiplier",
373 postSweepsLevelMultiplier_
374 );
375 controlDict_.readIfPresent("maxPostSweeps", maxPostSweeps_);
376 controlDict_.readIfPresent("nFinestSweeps", nFinestSweeps_);
377 controlDict_.readIfPresent("interpolateCorrection", interpolateCorrection_);
378 controlDict_.readIfPresent("scaleCorrection", scaleCorrection_);
379 controlDict_.readIfPresent("directSolveCoarsest", directSolveCoarsest_);
380
381 if ((log_ >= 2) || debug)
382 {
383 Info<< "GAMGSolver settings :"
384 << " cacheAgglomeration:" << cacheAgglomeration_
385 << " nPreSweeps:" << nPreSweeps_
386 << " preSweepsLevelMultiplier:" << preSweepsLevelMultiplier_
387 << " maxPreSweeps:" << maxPreSweeps_
388 << " nPostSweeps:" << nPostSweeps_
389 << " postSweepsLevelMultiplier:" << postSweepsLevelMultiplier_
390 << " maxPostSweeps:" << maxPostSweeps_
391 << " nFinestSweeps:" << nFinestSweeps_
392 << " interpolateCorrection:" << interpolateCorrection_
393 << " scaleCorrection:" << scaleCorrection_
394 << " directSolveCoarsest:" << directSolveCoarsest_
395 << endl;
396 }
397}
398
399
400const Foam::lduMatrix& Foam::GAMGSolver::matrixLevel(const label i) const
401{
402 return i ? matrixLevels_[i-1] : matrix_;
403}
404
405
406const Foam::lduInterfaceFieldPtrsList& Foam::GAMGSolver::interfaceLevel
407(
408 const label i
409) const
410{
411 return i ? interfaceLevels_[i-1] : interfaces_;
412}
413
414
415const Foam::FieldField<Foam::Field, Foam::scalar>&
416Foam::GAMGSolver::interfaceBouCoeffsLevel
417(
418 const label i
419) const
420{
421 return i ? interfaceLevelsBouCoeffs_[i-1] : interfaceBouCoeffs_;
422}
423
424
425const Foam::FieldField<Foam::Field, Foam::scalar>&
426Foam::GAMGSolver::interfaceIntCoeffsLevel
427(
428 const label i
429) const
430{
431 return i ? interfaceLevelsIntCoeffs_[i-1] : interfaceIntCoeffs_;
432}
433
434
435// ************************************************************************* //
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Geometric agglomerated algebraic multigrid agglomeration class.
bool processorAgglomerate() const
Whether to agglomerate across processors.
const labelList & procAgglomMap(const label fineLeveli) const
Mapping from processor to agglomerated processor (global, all processors have the same information)....
bool hasProcMesh(const label fineLeveli) const
Check that level has combined mesh.
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 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.
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.
Geometric agglomerated algebraic multigrid solver.
Definition GAMGSolver.H:75
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
virtual ~GAMGSolver()
Destructor.
Definition GAMGSolver.C:339
Similar to IStringStream but using an externally managed buffer for its input. This allows the input ...
void reset(const char *buffer, size_t nbytes)
Reset input area, position to buffer start and clear errors.
Class to perform the LU decomposition on a symmetric matrix.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An OSstream with internal List storage.
Preconditioned bi-conjugate gradient stabilized solver for asymmetric lduMatrices using a run-time se...
Definition PBiCGStab.H:66
Preconditioned conjugate gradient solver for symmetric lduMatrices using a run-time selectable precon...
Definition PCG.H:54
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
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:171
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
label size() const noexcept
The number of entries 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:133
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
virtual bool isDict() const noexcept
True if this entry is a dictionary.
Definition entry.H:284
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary, otherwise Fatal.
@ LITERAL_RECURSIVE
Definition keyType.H:87
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition lduMatrix.H:338
const lduMatrix & matrix_
Definition lduMatrix.H:158
solver(const word &fieldName, const lduMatrix &matrix, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, const dictionary &solverControls)
Construct solver for given field name, matrix etc.
scalar tolerance_
Final convergence tolerance.
Definition lduMatrix.H:191
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition lduMatrix.H:343
const lduMatrix & matrix() const noexcept
Definition lduMatrix.H:328
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 log_
Verbosity level for solver output statements.
Definition lduMatrix.H:171
scalar relTol_
Convergence tolerance relative to the initial.
Definition lduMatrix.H:196
virtual void readControls()
Read the control parameters from controlDict_.
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition lduMatrix.H:333
dictionary controlDict_
Dictionary of solution controls.
Definition lduMatrix.H:166
virtual const word & type() const =0
Runtime type information.
const word & fieldName() const noexcept
Definition lduMatrix.H:323
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
Base solver class.
Definition solver.H:48
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
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.
List< label > labelList
A List of labels.
Definition List.H:62
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
lduMatrix::solver::addasymMatrixConstructorToTable< GAMGSolver > addGAMGAsymSolverMatrixConstructorToTable_
Definition GAMGSolver.C:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
lduMatrix::solver::addsymMatrixConstructorToTable< GAMGSolver > addGAMGSolverMatrixConstructorToTable_
Definition GAMGSolver.C:34
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299