1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "GAMGSolver.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 defineTypeNameAndDebug(GAMGSolver, 0);
35 lduMatrix::solver::addsymMatrixConstructorToTable<GAMGSolver>
36 addGAMGSolverMatrixConstructorToTable_;
38 lduMatrix::solver::addasymMatrixConstructorToTable<GAMGSolver>
39 addGAMGAsymSolverMatrixConstructorToTable_;
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 Foam::GAMGSolver::GAMGSolver
47 const word& fieldName,
48 const lduMatrix& matrix,
49 const FieldField<Field, scalar>& interfaceBouCoeffs,
50 const FieldField<Field, scalar>& interfaceIntCoeffs,
51 const lduInterfaceFieldPtrsList& interfaces,
52 const dictionary& solverControls
65 // Default values for all controls
66 // which may be overridden by those in controlDict
67 cacheAgglomeration_(false),
71 scaleCorrection_(matrix.symmetric()),
72 directSolveCoarsest_(false),
73 agglomeration_(GAMGAgglomeration::New(matrix_, controlDict_)),
75 matrixLevels_(agglomeration_.size()),
76 interfaceLevels_(agglomeration_.size()),
77 interfaceLevelsBouCoeffs_(agglomeration_.size()),
78 interfaceLevelsIntCoeffs_(agglomeration_.size())
82 forAll(agglomeration_, fineLevelIndex)
84 agglomerateMatrix(fineLevelIndex);
87 if (matrixLevels_.size())
89 const label coarsestLevel = matrixLevels_.size() - 1;
91 if (directSolveCoarsest_)
93 coarsestLUMatrixPtr_.set
97 matrixLevels_[coarsestLevel],
98 interfaceLevelsBouCoeffs_[coarsestLevel],
99 interfaceLevels_[coarsestLevel]
108 "GAMGSolver::GAMGSolver"
110 "const word& fieldName,"
111 "const lduMatrix& matrix,"
112 "const FieldField<Field, scalar>& interfaceBouCoeffs,"
113 "const FieldField<Field, scalar>& interfaceIntCoeffs,"
114 "const lduInterfaceFieldPtrsList& interfaces,"
115 "const dictionary& solverControls"
117 ) << "No coarse levels created, either matrix too small for GAMG"
118 " or nCellsInCoarsestLevel too large.\n"
119 " Either choose another solver of reduce "
120 "nCellsInCoarsestLevel."
126 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
128 Foam::GAMGSolver::~GAMGSolver()
130 // Clear the the lists of pointers to the interfaces
131 forAll (interfaceLevels_, leveli)
133 lduInterfaceFieldPtrsList& curLevel = interfaceLevels_[leveli];
144 if (!cacheAgglomeration_)
146 delete &agglomeration_;
151 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
153 void Foam::GAMGSolver::readControls()
155 lduMatrix::solver::readControls();
157 // we could also consider supplying defaults here too
158 controlDict_.readIfPresent("cacheAgglomeration", cacheAgglomeration_);
159 controlDict_.readIfPresent("nPreSweeps", nPreSweeps_);
160 controlDict_.readIfPresent("nPostSweeps", nPostSweeps_);
161 controlDict_.readIfPresent("nFinestSweeps", nFinestSweeps_);
162 controlDict_.readIfPresent("scaleCorrection", scaleCorrection_);
163 controlDict_.readIfPresent("directSolveCoarsest", directSolveCoarsest_);
167 const Foam::lduMatrix& Foam::GAMGSolver::matrixLevel(const label i) const
175 return matrixLevels_[i - 1];
180 const Foam::lduInterfaceFieldPtrsList& Foam::GAMGSolver::interfaceLevel
191 return interfaceLevels_[i - 1];
196 const Foam::FieldField<Foam::Field, Foam::scalar>&
197 Foam::GAMGSolver::interfaceBouCoeffsLevel
204 return interfaceBouCoeffs_;
208 return interfaceLevelsBouCoeffs_[i - 1];
213 const Foam::FieldField<Foam::Field, Foam::scalar>&
214 Foam::GAMGSolver::interfaceIntCoeffsLevel
221 return interfaceIntCoeffs_;
225 return interfaceLevelsIntCoeffs_[i - 1];
230 // ************************************************************************* //