GAMGSolver.H
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-2016 OpenFOAM Foundation
9 Copyright (C) 2019 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
27Class
28 Foam::GAMGSolver
29
30Group
31 grpLduMatrixSolvers
32
33Description
34 Geometric agglomerated algebraic multigrid solver.
35
36 Characteristics:
37 - Requires positive definite, diagonally dominant matrix.
38 - Agglomeration algorithm: selectable and optionally cached.
39 - Restriction operator: summation.
40 - Prolongation operator: injection.
41 - Smoother: Gauss-Seidel.
42 - Coarse matrix creation: central coefficient: summation of fine grid
43 central coefficients with the removal of intra-cluster face;
44 off-diagonal coefficient: summation of off-diagonal faces.
45 - Coarse matrix scaling: performed by correction scaling, using steepest
46 descent optimisation.
47 - Type of cycle: V-cycle with optional pre-smoothing.
48 - Coarsest-level matrix solved using PCG or PBiCGStab.
49
50SourceFiles
51 GAMGSolver.C
52 GAMGSolverAgglomerateMatrix.C
53 GAMGSolverInterpolate.C
54 GAMGSolverScale.C
55 GAMGSolverSolve.C
56
57\*---------------------------------------------------------------------------*/
58
59#ifndef GAMGSolver_H
60#define GAMGSolver_H
61
62#include "GAMGAgglomeration.H"
63#include "lduMatrix.H"
64#include "labelField.H"
65#include "primitiveFields.H"
66#include "LUscalarMatrix.H"
67
68// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69
70namespace Foam
71{
72
73/*---------------------------------------------------------------------------*\
74 Class GAMGSolver Declaration
75\*---------------------------------------------------------------------------*/
77class GAMGSolver
78:
80{
81 // Private data
82
83 bool cacheAgglomeration_;
84
85 //- Number of pre-smoothing sweeps
86 label nPreSweeps_;
87
88 //- Lever multiplier for the number of pre-smoothing sweeps
89 label preSweepsLevelMultiplier_;
90
91 //- Maximum number of pre-smoothing sweeps
92 label maxPreSweeps_;
93
94 //- Number of post-smoothing sweeps
95 label nPostSweeps_;
96
97 //- Lever multiplier for the number of post-smoothing sweeps
98 label postSweepsLevelMultiplier_;
99
100 //- Maximum number of post-smoothing sweeps
101 label maxPostSweeps_;
102
103 //- Number of smoothing sweeps on finest mesh
104 label nFinestSweeps_;
105
106 //- Choose if the corrections should be interpolated after injection.
107 // By default corrections are not interpolated.
108 bool interpolateCorrection_;
109
110 //- Choose if the corrections should be scaled.
111 // By default corrections for symmetric matrices are scaled
112 // but not for asymmetric matrices.
113 bool scaleCorrection_;
114
115 //- Direct or iteratively solve the coarsest level
116 bool directSolveCoarsest_;
117
118 //- The agglomeration
119 const GAMGAgglomeration& agglomeration_;
120
121 //- Hierarchy of matrix levels
122 PtrList<lduMatrix> matrixLevels_;
123
124 //- Hierarchy of interfaces.
125 PtrList<PtrList<lduInterfaceField>> primitiveInterfaceLevels_;
126
127 //- Hierarchy of interfaces in lduInterfaceFieldPtrs form
128 PtrList<lduInterfaceFieldPtrsList> interfaceLevels_;
129
130 //- Hierarchy of interface boundary coefficients
131 PtrList<FieldField<Field, scalar>> interfaceLevelsBouCoeffs_;
132
133 //- Hierarchy of interface internal coefficients
134 PtrList<FieldField<Field, scalar>> interfaceLevelsIntCoeffs_;
135
136 //- LU decomposed coarsest matrix
137 autoPtr<LUscalarMatrix> coarsestLUMatrixPtr_;
138
139 //- Sparse coarsest matrix solver
140 autoPtr<lduMatrix::solver> coarsestSolverPtr_;
141
142
143 // Private Member Functions
144
145 //- Read control parameters from the control dictionary
146 virtual void readControls();
147
148 //- Simplified access to interface level
149 const lduInterfaceFieldPtrsList& interfaceLevel
150 (
151 const label i
152 ) const;
153
154 //- Simplified access to matrix level
155 const lduMatrix& matrixLevel(const label i) const;
156
157 //- Simplified access to interface boundary coeffs level
158 const FieldField<Field, scalar>& interfaceBouCoeffsLevel
159 (
160 const label i
161 ) const;
162
163 //- Simplified access to interface internal coeffs level
164 const FieldField<Field, scalar>& interfaceIntCoeffsLevel
165 (
166 const label i
167 ) const;
168
169 //- Agglomerate coarse matrix. Supply mesh to use - so we can
170 // construct temporary matrix on the fine mesh (instead of the coarse
171 // mesh)
172 void agglomerateMatrix
173 (
174 const label fineLevelIndex,
175 const lduMesh& coarseMesh,
176 const lduInterfacePtrsList& coarseMeshInterfaces
177 );
178
179 //- Agglomerate coarse interface coefficients
180 void agglomerateInterfaceCoefficients
181 (
182 const label fineLevelIndex,
183 const lduInterfacePtrsList& coarseMeshInterfaces,
184 PtrList<lduInterfaceField>& coarsePrimInterfaces,
185 lduInterfaceFieldPtrsList& coarseInterfaces,
186 FieldField<Field, scalar>& coarseInterfaceBouCoeffs,
187 FieldField<Field, scalar>& coarseInterfaceIntCoeffs
188 ) const;
189
190 //- Collect matrices from other processors
191 void gatherMatrices
192 (
193 const labelList& procIDs,
194 const lduMesh& dummyMesh,
195 const label meshComm,
196
197 const lduMatrix& mat,
201
202 PtrList<lduMatrix>& otherMats,
203 PtrList<FieldField<Field, scalar>>& otherBouCoeffs,
204 PtrList<FieldField<Field, scalar>>& otherIntCoeffs,
205 List<boolList>& otherTransforms,
206 List<List<label>>& otherRanks
207 ) const;
208
209 //- Agglomerate processor matrices
210 void procAgglomerateMatrix
211 (
212 // Agglomeration information
213 const labelList& procAgglomMap,
214 const List<label>& agglomProcIDs,
215
216 const label levelI,
217
218 // Resulting matrix
219 autoPtr<lduMatrix>& allMatrixPtr,
220 FieldField<Field, scalar>& allInterfaceBouCoeffs,
221 FieldField<Field, scalar>& allInterfaceIntCoeffs,
222 PtrList<lduInterfaceField>& allPrimitiveInterfaces,
223 lduInterfaceFieldPtrsList& allInterfaces
224 ) const;
225
226 //- Agglomerate processor matrices
227 void procAgglomerateMatrix
228 (
229 const labelList& procAgglomMap,
230 const List<label>& agglomProcIDs,
231 const label levelI
232 );
233
234 //- Interpolate the correction after injected prolongation
235 void interpolate
236 (
238 solveScalarField& Apsi,
239 const lduMatrix& m,
242 const direction cmpt
243 ) const;
244
245 //- Interpolate the correction after injected prolongation and
246 // re-normalise
247 void interpolate
248 (
250 solveScalarField& Apsi,
251 const lduMatrix& m,
254 const labelList& restrictAddressing,
255 const solveScalarField& psiC,
256 const direction cmpt
257 ) const;
258
259 //- Calculate and apply the scaling factor from Acf, coarseSource
260 // and coarseField.
261 // At the same time do a Jacobi iteration on the coarseField using
262 // the Acf provided after the coarseField values are used for the
263 // scaling factor.
264 void scale
265 (
267 solveScalarField& Acf,
268 const lduMatrix& A,
269 const FieldField<Field, scalar>& interfaceLevelBouCoeffs,
270 const lduInterfaceFieldPtrsList& interfaceLevel,
271 const solveScalarField& source,
272 const direction cmpt
273 ) const;
274
275 //- Initialise the data structures for the V-cycle
276 void initVcycle
277 (
278 PtrList<solveScalarField>& coarseCorrFields,
279 PtrList<solveScalarField>& coarseSources,
281 solveScalarField& scratch1,
282 solveScalarField& scratch2
283 ) const;
284
285
286 //- Perform a single GAMG V-cycle with pre, post and finest smoothing.
287 void Vcycle
288 (
289 const PtrList<lduMatrix::smoother>& smoothers,
291 const scalarField& source,
292 solveScalarField& Apsi,
293 solveScalarField& finestCorrection,
294 solveScalarField& finestResidual,
295
296 solveScalarField& scratch1,
297 solveScalarField& scratch2,
298
299 PtrList<solveScalarField>& coarseCorrFields,
300 PtrList<solveScalarField>& coarseSources,
301 const direction cmpt=0
302 ) const;
303
304 //- Create and return the dictionary to specify the PCG solver
305 // to solve the coarsest level
306 dictionary PCGsolverDict
307 (
308 const scalar tol,
309 const scalar relTol
310 ) const;
311
312 //- Create and return the dictionary to specify the PBiCGStab solver
313 // to solve the coarsest level
314 dictionary PBiCGStabSolverDict
315 (
316 const scalar tol,
317 const scalar relTol
318 ) const;
319
320 //- Solve the coarsest level with either an iterative or direct solver
321 void solveCoarsestLevel
322 (
323 solveScalarField& coarsestCorrField,
324 const solveScalarField& coarsestSource
325 ) const;
326
327
328public:
330 friend class GAMGPreconditioner;
331
332 //- Runtime type information
333 TypeName("GAMG");
334
335
336 // Constructors
337
338 //- Construct from lduMatrix and solver controls
340 (
341 const word& fieldName,
342 const lduMatrix& matrix,
346 const dictionary& solverControls
347 );
348
349
350 //- Destructor
351 virtual ~GAMGSolver();
352
353
354 // Member Functions
355
356 //- Solve
358 (
360 const scalarField& source,
361 const direction cmpt=0
362 ) const;
363};
364
365
366// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
367
368} // End namespace Foam
369
370// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
371
372#endif
373
374// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
Geometric agglomerated algebraic multigrid agglomeration class.
Geometric agglomerated algebraic multigrid preconditioner.
Geometric agglomerated algebraic multigrid solver.
Definition: GAMGSolver.H:79
TypeName("GAMG")
Runtime type information.
virtual ~GAMGSolver()
Destructor.
Definition: GAMGSolver.C:336
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:99
const FieldField< Field, scalar > & interfaceIntCoeffs() const noexcept
Definition: lduMatrix.H:243
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:248
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:233
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition: lduMatrix.H:238
const word & fieldName() const noexcept
Definition: lduMatrix.H:228
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:84
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:63
A class for handling words, derived from Foam::string.
Definition: word.H:68
const volScalarField & psi
rDeltaTY field()
Namespace for OpenFOAM.
uint8_t direction
Definition: direction.H:56
Specialisations of Field<T> for scalar, vector and tensor.
CEqn solve()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73