Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / OpenFOAM / matrices / lduMatrix / solvers / GAMG / GAMGSolver.H
blob41a660f3fcb4e50a219e69a8037bc66c1c79ab3a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Class
25     Foam::GAMGSolver
27 Description
28     Geometric agglomerated algebraic multigrid solver.
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.
44 SourceFiles
45     GAMGSolver.C
46     GAMGSolverCalcAgglomeration.C
47     GAMGSolverMakeCoarseMatrix.C
48     GAMGSolverOperations.C
49     GAMGSolverSolve.C
51 \*---------------------------------------------------------------------------*/
53 #ifndef GAMGSolver_H
54 #define GAMGSolver_H
56 #include "GAMGAgglomeration.H"
57 #include "lduMatrix.H"
58 #include "labelField.H"
59 #include "primitiveFields.H"
60 #include "LUscalarMatrix.H"
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 namespace Foam
67 /*---------------------------------------------------------------------------*\
68                            Class GAMGSolver Declaration
69 \*---------------------------------------------------------------------------*/
71 class GAMGSolver
73     public lduMatrix::solver
75     // Private data
77         bool cacheAgglomeration_;
79         //- Number of pre-smoothing sweeps
80         label nPreSweeps_;
82         //- Number of post-smoothing sweeps
83         label nPostSweeps_;
85         //- Number of smoothing sweeps on finest mesh
86         label nFinestSweeps_;
88         //- Choose if the corrections should be scaled.
89         //  By default corrections for symmetric matrices are scaled
90         //  but not for asymmetric matrices.
91         bool scaleCorrection_;
93         //- Direct or iteratively solve the coarsest level
94         bool directSolveCoarsest_;
96         //- The agglomeration
97         const GAMGAgglomeration& agglomeration_;
99         //- Hierarchy of matrix levels
100         PtrList<lduMatrix> matrixLevels_;
102         //- Hierarchy of interfaces.
103         //  Warning: Needs to be deleted explicitly.
104         PtrList<lduInterfaceFieldPtrsList> interfaceLevels_;
106         //- Hierarchy of interface boundary coefficients
107         PtrList<FieldField<Field, scalar> > interfaceLevelsBouCoeffs_;
109         //- Hierarchy of interface internal coefficients
110         PtrList<FieldField<Field, scalar> > interfaceLevelsIntCoeffs_;
112         //- LU decompsed coarsest matrix
113         autoPtr<LUscalarMatrix> coarsestLUMatrixPtr_;
116     // Private Member Functions
118         //- Read control parameters from the control dictionary
119         virtual void readControls();
121         //- Simplified access to interface level
122         const lduInterfaceFieldPtrsList& interfaceLevel
123         (
124             const label i
125         ) const;
127         //- Simplified access to matrix level
128         const lduMatrix& matrixLevel(const label i) const;
130         //- Simplified access to interface boundary coeffs level
131         const FieldField<Field, scalar>& interfaceBouCoeffsLevel
132         (
133             const label i
134         ) const;
136         //- Simplified access to interface internal coeffs level
137         const FieldField<Field, scalar>& interfaceIntCoeffsLevel
138         (
139             const label i
140         ) const;
142         //- Agglomerate coarse matrix
143         void agglomerateMatrix(const label fineLevelIndex);
145         //- Calculate and return the scaling factor from Acf, coarseSource
146         //  and coarseField.
147         //  At the same time do a Jacobi iteration on the coarseField using
148         //  the Acf provided after the coarseField values are used for the
149         //  scaling factor.
150         scalar scalingFactor
151         (
152             scalarField& field,
153             const scalarField& source,
154             const scalarField& Acf,
155             const scalarField& D
156         ) const;
158         //- Calculate Acf and calculate and return the scaling factor.
159         scalar scalingFactor
160         (
161             scalarField& Acf,
162             const lduMatrix& A,
163             scalarField& field,
164             const FieldField<Field, scalar>& interfaceLevelBouCoeffs,
165             const lduInterfaceFieldPtrsList& interfaceLevel,
166             const scalarField& source,
167             const direction cmpt
168         ) const;
171         //- Initialise the data structures for the V-cycle
172         void initVcycle
173         (
174             PtrList<scalarField>& coarseCorrFields,
175             PtrList<scalarField>& coarseSources,
176             PtrList<lduMatrix::smoother>& smoothers
177         ) const;
180         //- Perform a single GAMG V-cycle with pre, post and finest smoothing.
181         void Vcycle
182         (
183             const PtrList<lduMatrix::smoother>& smoothers,
184             scalarField& psi,
185             const scalarField& source,
186             scalarField& Apsi,
187             scalarField& finestCorrection,
188             scalarField& finestResidual,
189             PtrList<scalarField>& coarseCorrFields,
190             PtrList<scalarField>& coarseSources,
191             const direction cmpt=0
192         ) const;
195         //- Solve the coarsest level with either an iterative or direct solver
196         void solveCoarsestLevel
197         (
198             scalarField& coarsestCorrField,
199             const scalarField& coarsestSource
200         ) const;
203 public:
205     friend class GAMGPreconditioner;
207     //- Runtime type information
208     TypeName("GAMG");
211     // Constructors
213         //- Construct from lduMatrix and solver controls
214         GAMGSolver
215         (
216             const word& fieldName,
217             const lduMatrix& matrix,
218             const FieldField<Field, scalar>& interfaceBouCoeffs,
219             const FieldField<Field, scalar>& interfaceIntCoeffs,
220             const lduInterfaceFieldPtrsList& interfaces,
221             const dictionary& solverControls
222         );
225     //- Destructor
226     virtual ~GAMGSolver();
229     // Member Functions
231         //- Solve
232         virtual lduMatrix::solverPerformance solve
233         (
234             scalarField& psi,
235             const scalarField& source,
236             const direction cmpt=0
237         ) const;
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 } // End namespace Foam
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 #endif
249 // ************************************************************************* //