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 | Copyright (C) 2011-2014 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::GAMGSolver
26 
27 Description
28  Geometric agglomerated algebraic multigrid solver.
29 
30  Characteristics:
31  - Requires positive definite, diagonally dominant matrix.
32  - Agglomeration algorithm: selectable and optionally cached.
33  - Restriction operator: summation.
34  - Prolongation operator: injection.
35  - Smoother: Gauss-Seidel.
36  - Coarse matrix creation: central coefficient: summation of fine grid
37  central coefficients with the removal of intra-cluster face;
38  off-diagonal coefficient: summation of off-diagonal faces.
39  - Coarse matrix scaling: performed by correction scaling, using steepest
40  descent optimisation.
41  - Type of cycle: V-cycle with optional pre-smoothing.
42  - Coarsest-level matrix solved using ICCG or BICCG.
43 
44 SourceFiles
45  GAMGSolver.C
46  GAMGSolverAgglomerateMatrix.C
47  GAMGSolverInterpolate.C
48  GAMGSolverScale.C
49  GAMGSolverSolve.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef GAMGSolver_H
54 #define GAMGSolver_H
55 
56 #include "GAMGAgglomeration.H"
57 #include "lduMatrix.H"
58 #include "labelField.H"
59 #include "primitiveFields.H"
60 #include "LUscalarMatrix.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 namespace Foam
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class GAMGSolver Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 class GAMGSolver
72 :
73  public lduMatrix::solver
74 {
75  // Private data
76 
78 
79  //- Number of pre-smoothing sweeps
81 
82  //- Lever multiplier for the number of pre-smoothing sweeps
84 
85  //- Maximum number of pre-smoothing sweeps
87 
88  //- Number of post-smoothing sweeps
90 
91  //- Lever multiplier for the number of post-smoothing sweeps
93 
94  //- Maximum number of post-smoothing sweeps
96 
97  //- Number of smoothing sweeps on finest mesh
99 
100  //- Choose if the corrections should be interpolated after injection.
101  // By default corrections are not interpolated.
103 
104  //- Choose if the corrections should be scaled.
105  // By default corrections for symmetric matrices are scaled
106  // but not for asymmetric matrices.
107  bool scaleCorrection_;
108 
109  //- Direct or iteratively solve the coarsest level
111 
112  //- The agglomeration
114 
115  //- Hierarchy of matrix levels
117 
118  //- Hierarchy of interfaces.
120 
121  //- Hierarchy of interfaces in lduInterfaceFieldPtrs form
123 
124  //- Hierarchy of interface boundary coefficients
126 
127  //- Hierarchy of interface internal coefficients
129 
130  //- LU decompsed coarsest matrix
132 
133 
134  // Private Member Functions
135 
136  //- Read control parameters from the control dictionary
137  virtual void readControls();
138 
139  //- Simplified access to interface level
141  (
142  const label i
143  ) const;
144 
145  //- Simplified access to matrix level
146  const lduMatrix& matrixLevel(const label i) const;
147 
148  //- Simplified access to interface boundary coeffs level
150  (
151  const label i
152  ) const;
153 
154  //- Simplified access to interface internal coeffs level
156  (
157  const label i
158  ) const;
159 
160  //- Agglomerate coarse matrix. Supply mesh to use - so we can
161  // construct temporary matrix on the fine mesh (instead of the coarse
162  // mesh)
163  void agglomerateMatrix
164  (
165  const label fineLevelIndex,
166  const lduMesh& coarseMesh,
167  const lduInterfacePtrsList& coarseMeshInterfaces
168  );
169 
170  //- Agglomerate coarse interface coefficients
172  (
173  const label fineLevelIndex,
174  const lduInterfacePtrsList& coarseMeshInterfaces,
175  PtrList<lduInterfaceField>& coarsePrimInterfaces,
176  lduInterfaceFieldPtrsList& coarseInterfaces,
177  FieldField<Field, scalar>& coarseInterfaceBouCoeffs,
178  FieldField<Field, scalar>& coarseInterfaceIntCoeffs
179  ) const;
180 
181  //- Collect matrices from other processors
182  void gatherMatrices
183  (
184  const labelList& procIDs,
185  const lduMesh& dummyMesh,
186  const label meshComm,
187 
188  const lduMatrix& mat,
192 
193  PtrList<lduMatrix>& otherMats,
194  PtrList<FieldField<Field, scalar> >& otherBouCoeffs,
195  PtrList<FieldField<Field, scalar> >& otherIntCoeffs,
196  List<boolList>& otherTransforms,
197  List<List<label> >& otherRanks
198  ) const;
199 
200  //- Agglomerate processor matrices
202  (
203  // Agglomeration information
204  const labelList& procAgglomMap,
205  const List<label>& agglomProcIDs,
206 
207  const label levelI,
208 
209  // Resulting matrix
210  autoPtr<lduMatrix>& allMatrixPtr,
211  FieldField<Field, scalar>& allInterfaceBouCoeffs,
212  FieldField<Field, scalar>& allInterfaceIntCoeffs,
213  PtrList<lduInterfaceField>& allPrimitiveInterfaces,
214  lduInterfaceFieldPtrsList& allInterfaces
215  ) const;
216 
217  //- Agglomerate processor matrices
219  (
220  const labelList& procAgglomMap,
221  const List<label>& agglomProcIDs,
222  const label levelI
223  );
224 
225  //- Interpolate the correction after injected prolongation
226  void interpolate
227  (
228  scalarField& psi,
229  scalarField& Apsi,
230  const lduMatrix& m,
233  const direction cmpt
234  ) const;
235 
236  //- Interpolate the correction after injected prolongation and
238  void interpolate
239  (
240  scalarField& psi,
241  scalarField& Apsi,
242  const lduMatrix& m,
245  const labelList& restrictAddressing,
246  const scalarField& psiC,
247  const direction cmpt
248  ) const;
249 
250  //- Calculate and apply the scaling factor from Acf, coarseSource
251  // and coarseField.
252  // At the same time do a Jacobi iteration on the coarseField using
253  // the Acf provided after the coarseField values are used for the
254  // scaling factor.
255  void scale
256  (
257  scalarField& field,
258  scalarField& Acf,
259  const lduMatrix& A,
260  const FieldField<Field, scalar>& interfaceLevelBouCoeffs,
262  const scalarField& source,
263  const direction cmpt
264  ) const;
265 
266  //- Initialise the data structures for the V-cycle
267  void initVcycle
268  (
269  PtrList<scalarField>& coarseCorrFields,
270  PtrList<scalarField>& coarseSources,
271  PtrList<lduMatrix::smoother>& smoothers,
272  scalarField& scratch1,
273  scalarField& scratch2
274  ) const;
275 
276 
277  //- Perform a single GAMG V-cycle with pre, post and finest smoothing.
278  void Vcycle
279  (
280  const PtrList<lduMatrix::smoother>& smoothers,
281  scalarField& psi,
282  const scalarField& source,
283  scalarField& Apsi,
284  scalarField& finestCorrection,
285  scalarField& finestResidual,
286 
287  scalarField& scratch1,
288  scalarField& scratch2,
289 
290  PtrList<scalarField>& coarseCorrFields,
291  PtrList<scalarField>& coarseSources,
292  const direction cmpt=0
293  ) const;
294 
295 
296  //- Solve the coarsest level with either an iterative or direct solver
297  void solveCoarsestLevel
298  (
299  scalarField& coarsestCorrField,
300  const scalarField& coarsestSource
301  ) const;
302 
303 
304 public:
305 
306  friend class GAMGPreconditioner;
307 
308  //- Runtime type information
309  TypeName("GAMG");
310 
311 
312  // Constructors
313 
314  //- Construct from lduMatrix and solver controls
315  GAMGSolver
316  (
317  const word& fieldName,
318  const lduMatrix& matrix,
322  const dictionary& solverControls
323  );
324 
325 
326  //- Destructor
327  virtual ~GAMGSolver();
328 
329 
330  // Member Functions
331 
332  //- Solve
333  virtual solverPerformance solve
334  (
335  scalarField& psi,
336  const scalarField& source,
337  const direction cmpt=0
338  ) const;
339 };
340 
341 
342 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
343 
344 } // End namespace Foam
345 
346 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347 
348 #endif
349 
350 // ************************************************************************* //
Foam::GAMGSolver::postSweepsLevelMultiplier_
label postSweepsLevelMultiplier_
Lever multiplier for the number of post-smoothing sweeps.
Definition: GAMGSolver.H:91
Foam::lduMatrix::solver::interfaces
const lduInterfaceFieldPtrsList & interfaces() const
Definition: lduMatrix.H:236
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
Foam::GAMGAgglomeration
Geometric agglomerated algebraic multigrid agglomeration class.
Definition: GAMGAgglomeration.H:61
Foam::GAMGSolver::initVcycle
void initVcycle(PtrList< scalarField > &coarseCorrFields, PtrList< scalarField > &coarseSources, PtrList< lduMatrix::smoother > &smoothers, scalarField &scratch1, scalarField &scratch2) const
Initialise the data structures for the V-cycle.
Definition: GAMGSolverSolve.C:449
Foam::GAMGSolver::preSweepsLevelMultiplier_
label preSweepsLevelMultiplier_
Lever multiplier for the number of pre-smoothing sweeps.
Definition: GAMGSolver.H:82
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:77
Foam::GAMGSolver::agglomerateInterfaceCoefficients
void agglomerateInterfaceCoefficients(const label fineLevelIndex, const lduInterfacePtrsList &coarseMeshInterfaces, PtrList< lduInterfaceField > &coarsePrimInterfaces, lduInterfaceFieldPtrsList &coarseInterfaces, FieldField< Field, scalar > &coarseInterfaceBouCoeffs, FieldField< Field, scalar > &coarseInterfaceIntCoeffs) const
Agglomerate coarse interface coefficients.
Definition: GAMGSolverAgglomerateMatrix.C:198
Foam::GAMGSolver::coarsestLUMatrixPtr_
autoPtr< LUscalarMatrix > coarsestLUMatrixPtr_
LU decompsed coarsest matrix.
Definition: GAMGSolver.H:130
Foam::GAMGSolver::directSolveCoarsest_
bool directSolveCoarsest_
Direct or iteratively solve the coarsest level.
Definition: GAMGSolver.H:109
Foam::GAMGSolver::solve
virtual solverPerformance solve(scalarField &psi, const scalarField &source, const direction cmpt=0) const
Solve.
Definition: GAMGSolverSolve.C:34
Foam::GAMGSolver::matrixLevel
const lduMatrix & matrixLevel(const label i) const
Simplified access to matrix level.
Definition: GAMGSolver.C:345
Foam::lduMatrix::solver::matrix
const lduMatrix & matrix() const
Definition: lduMatrix.H:221
Foam::GAMGSolver::interfaceLevelsIntCoeffs_
PtrList< FieldField< Field, scalar > > interfaceLevelsIntCoeffs_
Hierarchy of interface internal coefficients.
Definition: GAMGSolver.H:127
Foam::GAMGSolver::~GAMGSolver
virtual ~GAMGSolver()
Destructor.
Definition: GAMGSolver.C:291
GAMGAgglomeration.H
Foam::lduMatrix::solver
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:91
lduMatrix.H
LUscalarMatrix.H
primitiveFields.H
Specialisations of Field<T> for scalar, vector and tensor.
Foam::GAMGSolver::interfaceLevel
const lduInterfaceFieldPtrsList & interfaceLevel(const label i) const
Simplified access to interface level.
Definition: GAMGSolver.C:359
A
simpleMatrix< scalar > A(Nc)
Foam::GAMGSolver::readControls
virtual void readControls()
Read control parameters from the control dictionary.
Definition: GAMGSolver.C:302
Foam::GAMGSolver::Vcycle
void Vcycle(const PtrList< lduMatrix::smoother > &smoothers, scalarField &psi, const scalarField &source, scalarField &Apsi, scalarField &finestCorrection, scalarField &finestResidual, scalarField &scratch1, scalarField &scratch2, PtrList< scalarField > &coarseCorrFields, PtrList< scalarField > &coarseSources, const direction cmpt=0) const
Perform a single GAMG V-cycle with pre, post and finest smoothing.
Definition: GAMGSolverSolve.C:151
Foam::GAMGSolver::matrixLevels_
PtrList< lduMatrix > matrixLevels_
Hierarchy of matrix levels.
Definition: GAMGSolver.H:115
Foam::GAMGSolver::interfaceLevelsBouCoeffs_
PtrList< FieldField< Field, scalar > > interfaceLevelsBouCoeffs_
Hierarchy of interface boundary coefficients.
Definition: GAMGSolver.H:124
Foam::GAMGSolver::interfaceIntCoeffsLevel
const FieldField< Field, scalar > & interfaceIntCoeffsLevel(const label i) const
Simplified access to interface internal coeffs level.
Definition: GAMGSolver.C:393
Foam::GAMGSolver::maxPostSweeps_
label maxPostSweeps_
Maximum number of post-smoothing sweeps.
Definition: GAMGSolver.H:94
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::GAMGSolver::nPostSweeps_
label nPostSweeps_
Number of post-smoothing sweeps.
Definition: GAMGSolver.H:88
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::GAMGSolver::scaleCorrection_
bool scaleCorrection_
Choose if the corrections should be scaled.
Definition: GAMGSolver.H:106
Foam::GAMGSolver::nPreSweeps_
label nPreSweeps_
Number of pre-smoothing sweeps.
Definition: GAMGSolver.H:79
Foam::UPtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:53
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::GAMGSolver::interpolate
void interpolate(scalarField &psi, scalarField &Apsi, const lduMatrix &m, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
Interpolate the correction after injected prolongation.
Definition: GAMGSolverInterpolate.C:31
Foam::GAMGSolver::cacheAgglomeration_
bool cacheAgglomeration_
Definition: GAMGSolver.H:76
Foam::GAMGSolver::TypeName
TypeName("GAMG")
Runtime type information.
Foam::GAMGSolver::solveCoarsestLevel
void solveCoarsestLevel(scalarField &coarsestCorrField, const scalarField &coarsestSource) const
Solve the coarsest level with either an iterative or direct solver.
Definition: GAMGSolverSolve.C:523
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::GAMGSolver::interfaceLevels_
PtrList< lduInterfaceFieldPtrsList > interfaceLevels_
Hierarchy of interfaces in lduInterfaceFieldPtrs form.
Definition: GAMGSolver.H:121
Foam::GAMGSolver::primitiveInterfaceLevels_
PtrList< PtrList< lduInterfaceField > > primitiveInterfaceLevels_
Hierarchy of interfaces.
Definition: GAMGSolver.H:118
Foam::GAMGPreconditioner
Geometric agglomerated algebraic multigrid preconditioner.
Definition: GAMGPreconditioner.H:51
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::GAMGSolver::interpolateCorrection_
bool interpolateCorrection_
Choose if the corrections should be interpolated after injection.
Definition: GAMGSolver.H:101
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::GAMGSolver::procAgglomerateMatrix
void procAgglomerateMatrix(const labelList &procAgglomMap, const List< label > &agglomProcIDs, const label levelI, autoPtr< lduMatrix > &allMatrixPtr, FieldField< Field, scalar > &allInterfaceBouCoeffs, FieldField< Field, scalar > &allInterfaceIntCoeffs, PtrList< lduInterfaceField > &allPrimitiveInterfaces, lduInterfaceFieldPtrsList &allInterfaces) const
Agglomerate processor matrices.
Definition: GAMGSolverAgglomerateMatrix.C:415
Foam::GAMGSolver::maxPreSweeps_
label maxPreSweeps_
Maximum number of pre-smoothing sweeps.
Definition: GAMGSolver.H:85
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::GAMGSolver
Geometric agglomerated algebraic multigrid solver.
Definition: GAMGSolver.H:70
Foam::GAMGSolver::agglomeration_
const GAMGAgglomeration & agglomeration_
The agglomeration.
Definition: GAMGSolver.H:112
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::GAMGSolver::interfaceBouCoeffsLevel
const FieldField< Field, scalar > & interfaceBouCoeffsLevel(const label i) const
Simplified access to interface boundary coeffs level.
Definition: GAMGSolver.C:376
Foam::lduMatrix::solver::interfaceBouCoeffs
const FieldField< Field, scalar > & interfaceBouCoeffs() const
Definition: lduMatrix.H:226
Foam::GAMGSolver::scale
void scale(scalarField &field, scalarField &Acf, const lduMatrix &A, const FieldField< Field, scalar > &interfaceLevelBouCoeffs, const lduInterfaceFieldPtrsList &interfaceLevel, const scalarField &source, const direction cmpt) const
Calculate and apply the scaling factor from Acf, coarseSource.
Definition: GAMGSolverScale.C:32
Foam::GAMGSolver::GAMGSolver
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:46
Foam::GAMGSolver::gatherMatrices
void gatherMatrices(const labelList &procIDs, const lduMesh &dummyMesh, const label meshComm, const lduMatrix &mat, const FieldField< Field, scalar > &interfaceBouCoeffs, const FieldField< Field, scalar > &interfaceIntCoeffs, const lduInterfaceFieldPtrsList &interfaces, PtrList< lduMatrix > &otherMats, PtrList< FieldField< Field, scalar > > &otherBouCoeffs, PtrList< FieldField< Field, scalar > > &otherIntCoeffs, List< boolList > &otherTransforms, List< List< label > > &otherRanks) const
Collect matrices from other processors.
Definition: GAMGSolverAgglomerateMatrix.C:285
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::lduMatrix::solver::fieldName
const word & fieldName() const
Definition: lduMatrix.H:216
Foam::lduMatrix::solver::interfaceIntCoeffs
const FieldField< Field, scalar > & interfaceIntCoeffs() const
Definition: lduMatrix.H:231
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:49
Foam::GAMGSolver::agglomerateMatrix
void agglomerateMatrix(const label fineLevelIndex, const lduMesh &coarseMesh, const lduInterfacePtrsList &coarseMeshInterfaces)
Agglomerate coarse matrix. Supply mesh to use - so we can.
Definition: GAMGSolverAgglomerateMatrix.C:34
labelField.H
Foam::lduMesh
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:51
Foam::GAMGSolver::nFinestSweeps_
label nFinestSweeps_
Number of smoothing sweeps on finest mesh.
Definition: GAMGSolver.H:97