GAMGSolver.C
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-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "GAMGSolver.H"
29 #include "GAMGInterface.H"
30 #include "PCG.H"
31 #include "PBiCGStab.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(GAMGSolver, 0);
38 
39  lduMatrix::solver::addsymMatrixConstructorToTable<GAMGSolver>
41 
42  lduMatrix::solver::addasymMatrixConstructorToTable<GAMGSolver>
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& fieldName,
52  const lduMatrix& matrix,
53  const FieldField<Field, scalar>& interfaceBouCoeffs,
54  const FieldField<Field, scalar>& interfaceIntCoeffs,
55  const lduInterfaceFieldPtrsList& interfaces,
56  const dictionary& solverControls
57 )
58 :
60  (
61  fieldName,
62  matrix,
63  interfaceBouCoeffs,
64  interfaceIntCoeffs,
65  interfaces,
66  solverControls
67  ),
68 
69  // Default values for all controls
70  // which may be overridden by those in controlDict
71  cacheAgglomeration_(true),
72  nPreSweeps_(0),
73  preSweepsLevelMultiplier_(1),
74  maxPreSweeps_(4),
75  nPostSweeps_(2),
76  postSweepsLevelMultiplier_(1),
77  maxPostSweeps_(4),
78  nFinestSweeps_(2),
79  interpolateCorrection_(false),
80  scaleCorrection_(matrix.symmetric()),
81  directSolveCoarsest_(false),
82  agglomeration_(GAMGAgglomeration::New(matrix_, controlDict_)),
83 
84  matrixLevels_(agglomeration_.size()),
85  primitiveInterfaceLevels_(agglomeration_.size()),
86  interfaceLevels_(agglomeration_.size()),
87  interfaceLevelsBouCoeffs_(agglomeration_.size()),
88  interfaceLevelsIntCoeffs_(agglomeration_.size())
89 {
90  readControls();
91 
92  if (agglomeration_.processorAgglomerate())
93  {
94  forAll(agglomeration_, fineLevelIndex)
95  {
96  if (agglomeration_.hasMeshLevel(fineLevelIndex))
97  {
98  if
99  (
100  (fineLevelIndex+1) < agglomeration_.size()
101  && agglomeration_.hasProcMesh(fineLevelIndex+1)
102  )
103  {
104  // Construct matrix without referencing the coarse mesh so
105  // construct a dummy mesh instead. This will get overwritten
106  // by the call to procAgglomerateMatrix so is only to get
107  // it through agglomerateMatrix
108 
109 
110  const lduInterfacePtrsList& fineMeshInterfaces =
111  agglomeration_.interfaceLevel(fineLevelIndex);
112 
113  PtrList<GAMGInterface> dummyPrimMeshInterfaces
114  (
115  fineMeshInterfaces.size()
116  );
117  lduInterfacePtrsList dummyMeshInterfaces
118  (
119  dummyPrimMeshInterfaces.size()
120  );
121  forAll(fineMeshInterfaces, intI)
122  {
123  if (fineMeshInterfaces.set(intI))
124  {
125  OStringStream os;
126  refCast<const GAMGInterface>
127  (
128  fineMeshInterfaces[intI]
129  ).write(os);
130  IStringStream is(os.str());
131 
132  dummyPrimMeshInterfaces.set
133  (
134  intI,
136  (
137  fineMeshInterfaces[intI].type(),
138  intI,
139  dummyMeshInterfaces,
140  is
141  )
142  );
143  }
144  }
145 
146  forAll(dummyPrimMeshInterfaces, intI)
147  {
148  if (dummyPrimMeshInterfaces.set(intI))
149  {
150  dummyMeshInterfaces.set
151  (
152  intI,
153  &dummyPrimMeshInterfaces[intI]
154  );
155  }
156  }
157 
158  // So:
159  // - pass in incorrect mesh (= fine mesh instead of coarse)
160  // - pass in dummy interfaces
161  agglomerateMatrix
162  (
163  fineLevelIndex,
164  agglomeration_.meshLevel(fineLevelIndex),
165  dummyMeshInterfaces
166  );
167 
168 
169  const labelList& procAgglomMap =
170  agglomeration_.procAgglomMap(fineLevelIndex+1);
171  const List<label>& procIDs =
172  agglomeration_.agglomProcIDs(fineLevelIndex+1);
173 
174  procAgglomerateMatrix
175  (
176  procAgglomMap,
177  procIDs,
178  fineLevelIndex
179  );
180  }
181  else
182  {
183  agglomerateMatrix
184  (
185  fineLevelIndex,
186  agglomeration_.meshLevel(fineLevelIndex + 1),
187  agglomeration_.interfaceLevel(fineLevelIndex + 1)
188  );
189  }
190  }
191  else
192  {
193  // No mesh. Not involved in calculation anymore
194  }
195  }
196  }
197  else
198  {
199  forAll(agglomeration_, fineLevelIndex)
200  {
201  // Agglomerate on to coarse level mesh
202  agglomerateMatrix
203  (
204  fineLevelIndex,
205  agglomeration_.meshLevel(fineLevelIndex + 1),
206  agglomeration_.interfaceLevel(fineLevelIndex + 1)
207  );
208  }
209  }
210 
211 
212  if (debug & 2)
213  {
214  for
215  (
216  label fineLevelIndex = 0;
217  fineLevelIndex <= matrixLevels_.size();
218  fineLevelIndex++
219  )
220  {
221  if (fineLevelIndex == 0 || matrixLevels_.set(fineLevelIndex-1))
222  {
223  const lduMatrix& matrix = matrixLevel(fineLevelIndex);
224  const lduInterfaceFieldPtrsList& interfaces =
225  interfaceLevel(fineLevelIndex);
226 
227  Pout<< "level:" << fineLevelIndex << nl
228  << " nCells:" << matrix.diag().size() << nl
229  << " nFaces:" << matrix.lower().size() << nl
230  << " nInterfaces:" << interfaces.size()
231  << endl;
232 
233  forAll(interfaces, i)
234  {
235  if (interfaces.set(i))
236  {
237  Pout<< " " << i
238  << "\ttype:" << interfaces[i].type()
239  << endl;
240  }
241  }
242  }
243  else
244  {
245  Pout<< "level:" << fineLevelIndex << " : no matrix" << endl;
246  }
247  }
248  Pout<< endl;
249  }
250 
251 
252  if (matrixLevels_.size())
253  {
254  const label coarsestLevel = matrixLevels_.size() - 1;
255 
256  if (matrixLevels_.set(coarsestLevel))
257  {
258  if (directSolveCoarsest_)
259  {
260  coarsestLUMatrixPtr_.reset
261  (
262  new LUscalarMatrix
263  (
264  matrixLevels_[coarsestLevel],
265  interfaceLevelsBouCoeffs_[coarsestLevel],
266  interfaceLevels_[coarsestLevel]
267  )
268  );
269  }
270  else
271  {
272  entry* coarseEntry = controlDict_.findEntry
273  (
274  "coarsestLevelCorr",
275  keyType::LITERAL_RECURSIVE
276  );
277  if (coarseEntry && coarseEntry->isDict())
278  {
279  coarsestSolverPtr_ = lduMatrix::solver::New
280  (
281  "coarsestLevelCorr",
282  matrixLevels_[coarsestLevel],
283  interfaceLevelsBouCoeffs_[coarsestLevel],
284  interfaceLevelsIntCoeffs_[coarsestLevel],
285  interfaceLevels_[coarsestLevel],
286  coarseEntry->dict()
287  );
288  }
289  else if (matrixLevels_[coarsestLevel].asymmetric())
290  {
291  coarsestSolverPtr_.reset
292  (
293  new PBiCGStab
294  (
295  "coarsestLevelCorr",
296  matrixLevels_[coarsestLevel],
297  interfaceLevelsBouCoeffs_[coarsestLevel],
298  interfaceLevelsIntCoeffs_[coarsestLevel],
299  interfaceLevels_[coarsestLevel],
300  PBiCGStabSolverDict(tolerance_, relTol_)
301  )
302  );
303  }
304  else
305  {
306  coarsestSolverPtr_.reset
307  (
308  new PCG
309  (
310  "coarsestLevelCorr",
311  matrixLevels_[coarsestLevel],
312  interfaceLevelsBouCoeffs_[coarsestLevel],
313  interfaceLevelsIntCoeffs_[coarsestLevel],
314  interfaceLevels_[coarsestLevel],
315  PCGsolverDict(tolerance_, relTol_)
316  )
317  );
318  }
319  }
320  }
321  }
322  else
323  {
325  << "No coarse levels created, either matrix too small for GAMG"
326  " or nCellsInCoarsestLevel too large.\n"
327  " Either choose another solver of reduce "
328  "nCellsInCoarsestLevel."
329  << exit(FatalError);
330  }
331 }
332 
333 
334 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
335 
337 {
338  if (!cacheAgglomeration_)
339  {
340  delete &agglomeration_;
341  }
342 }
343 
344 
345 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
346 
347 void Foam::GAMGSolver::readControls()
348 {
350 
351  controlDict_.readIfPresent("cacheAgglomeration", cacheAgglomeration_);
352  controlDict_.readIfPresent("nPreSweeps", nPreSweeps_);
353  controlDict_.readIfPresent
354  (
355  "preSweepsLevelMultiplier",
356  preSweepsLevelMultiplier_
357  );
358  controlDict_.readIfPresent("maxPreSweeps", maxPreSweeps_);
359  controlDict_.readIfPresent("nPostSweeps", nPostSweeps_);
360  controlDict_.readIfPresent
361  (
362  "postSweepsLevelMultiplier",
363  postSweepsLevelMultiplier_
364  );
365  controlDict_.readIfPresent("maxPostSweeps", maxPostSweeps_);
366  controlDict_.readIfPresent("nFinestSweeps", nFinestSweeps_);
367  controlDict_.readIfPresent("interpolateCorrection", interpolateCorrection_);
368  controlDict_.readIfPresent("scaleCorrection", scaleCorrection_);
369  controlDict_.readIfPresent("directSolveCoarsest", directSolveCoarsest_);
370 
371  if (debug)
372  {
373  Info<< "GAMGSolver settings :"
374  << " cacheAgglomeration:" << cacheAgglomeration_
375  << " nPreSweeps:" << nPreSweeps_
376  << " preSweepsLevelMultiplier:" << preSweepsLevelMultiplier_
377  << " maxPreSweeps:" << maxPreSweeps_
378  << " nPostSweeps:" << nPostSweeps_
379  << " postSweepsLevelMultiplier:" << postSweepsLevelMultiplier_
380  << " maxPostSweeps:" << maxPostSweeps_
381  << " nFinestSweeps:" << nFinestSweeps_
382  << " interpolateCorrection:" << interpolateCorrection_
383  << " scaleCorrection:" << scaleCorrection_
384  << " directSolveCoarsest:" << directSolveCoarsest_
385  << endl;
386  }
387 }
388 
389 
390 const Foam::lduMatrix& Foam::GAMGSolver::matrixLevel(const label i) const
391 {
392  if (i == 0)
393  {
394  return matrix_;
395  }
396  else
397  {
398  return matrixLevels_[i - 1];
399  }
400 }
401 
402 
403 const Foam::lduInterfaceFieldPtrsList& Foam::GAMGSolver::interfaceLevel
404 (
405  const label i
406 ) const
407 {
408  if (i == 0)
409  {
410  return interfaces_;
411  }
412  else
413  {
414  return interfaceLevels_[i - 1];
415  }
416 }
417 
418 
420 Foam::GAMGSolver::interfaceBouCoeffsLevel
421 (
422  const label i
423 ) const
424 {
425  if (i == 0)
426  {
427  return interfaceBouCoeffs_;
428  }
429  else
430  {
431  return interfaceLevelsBouCoeffs_[i - 1];
432  }
433 }
434 
435 
437 Foam::GAMGSolver::interfaceIntCoeffsLevel
438 (
439  const label i
440 ) const
441 {
442  if (i == 0)
443  {
444  return interfaceIntCoeffs_;
445  }
446  else
447  {
448  return interfaceLevelsIntCoeffs_[i - 1];
449  }
450 }
451 
452 
453 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::UPtrList::size
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:97
Foam::OSstream::write
virtual bool write(const token &tok)
Write token to stream or otherwise handle it.
Definition: OSstream.C:36
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::addGAMGAsymSolverMatrixConstructorToTable_
lduMatrix::solver::addasymMatrixConstructorToTable< GAMGSolver > addGAMGAsymSolverMatrixConstructorToTable_
Definition: GAMGSolver.C:43
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::PCG
Preconditioned conjugate gradient solver for symmetric lduMatrices using a run-time selectable precon...
Definition: PCG.H:55
Foam::addGAMGSolverMatrixConstructorToTable_
lduMatrix::solver::addsymMatrixConstructorToTable< GAMGSolver > addGAMGSolverMatrixConstructorToTable_
Definition: GAMGSolver.C:40
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:83
Foam::lduMatrix::symmetric
bool symmetric() const
Definition: lduMatrix.H:622
Foam::GAMGSolver::~GAMGSolver
virtual ~GAMGSolver()
Destructor.
Definition: GAMGSolver.C:336
Foam::lduMatrix::solver
Abstract base-class for lduMatrix solvers.
Definition: lduMatrix.H:97
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
GAMGInterface.H
Foam::entry::isDict
virtual bool isDict() const
Return true if this entry is a dictionary.
Definition: entry.H:222
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::UPtrList< const lduInterfaceField >
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::lduMatrix::solver::readControls
virtual void readControls()
Read the control parameters from the controlDict_.
Definition: lduMatrixSolver.C:163
Foam::LUscalarMatrix
Class to perform the LU decomposition on a symmetric matrix.
Definition: LUscalarMatrix.H:56
Foam::IStringStream
Input from string buffer, using a ISstream.
Definition: StringStream.H:111
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::entry::dict
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary.
Foam::PBiCGStab
Preconditioned bi-conjugate gradient stabilized solver for asymmetric lduMatrices using a run-time se...
Definition: PBiCGStab.H:67
Foam::lduMatrix::diag
scalarField & diag()
Definition: lduMatrix.C:192
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::Detail::StringStreamAllocator::str
Foam::string str() const
Get the string - as Foam::string rather than std::string.
Definition: StringStream.H:91
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::UPtrList::set
const T * set(const label i) const
Return const pointer to element (can be nullptr),.
Definition: UPtrList.H:170
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::List< label >
Foam::OStringStream
Output to string buffer, using a OSstream.
Definition: StringStream.H:196
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:50
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::lduMatrix::lower
scalarField & lower()
Definition: lduMatrix.C:174
PCG.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
GAMGSolver.H
PBiCGStab.H