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