MGridGenGAMGAgglomeration.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  Copyright (C) 2018 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 \*---------------------------------------------------------------------------*/
28 
30 #include "fvMesh.H"
32 #include "processorLduInterface.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(MGridGenGAMGAgglomeration, 0);
39 
41  (
42  GAMGAgglomeration,
43  MGridGenGAMGAgglomeration,
44  lduMesh
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::MGridGenGAMGAgglomeration::swap
52 (
53  const lduInterfacePtrsList& interfaces,
54  const labelUList& cellValues,
55  PtrList<labelList>& nbrValues
56 ) const
57 {
58  const label nReq = Pstream::nRequests();
59 
60  // Initialise transfer of restrict addressing on the interface
61  forAll(interfaces, inti)
62  {
63  if (interfaces.set(inti))
64  {
65  interfaces[inti].initInternalFieldTransfer
66  (
68  cellValues
69  );
70  }
71  }
72 
73  if (Pstream::parRun())
74  {
76  }
77 
78  // Get the interface agglomeration
79  nbrValues.setSize(interfaces.size());
80  forAll(interfaces, inti)
81  {
82  if (interfaces.set(inti))
83  {
84  nbrValues.set
85  (
86  inti,
87  new labelList
88  (
89  interfaces[inti].internalFieldTransfer
90  (
92  cellValues
93  )
94  )
95  );
96  }
97  }
98 }
99 
100 
101 void Foam::MGridGenGAMGAgglomeration::getNbrAgglom
102 (
103  const lduAddressing& addr,
104  const lduInterfacePtrsList& interfaces,
105  const PtrList<labelList>& nbrGlobalAgglom,
106  labelList& cellToNbrAgglom
107 ) const
108 {
109  cellToNbrAgglom.setSize(addr.size());
110  cellToNbrAgglom = -1;
111 
112  forAll(interfaces, inti)
113  {
114  if (interfaces.set(inti))
115  {
116  if (isA<processorLduInterface>(interfaces[inti]))
117  {
118  const processorLduInterface& pldui =
119  refCast<const processorLduInterface>(interfaces[inti]);
120 
121  if (pldui.myProcNo() > pldui.neighbProcNo())
122  {
123  const labelUList& faceCells =
124  interfaces[inti].faceCells();
125  const labelList& nbrData = nbrGlobalAgglom[inti];
126 
127  forAll(faceCells, i)
128  {
129  cellToNbrAgglom[faceCells[i]] = nbrData[i];
130  }
131  }
132  }
133  }
134  }
135 }
136 
137 
138 void Foam::MGridGenGAMGAgglomeration::detectSharedFaces
139 (
140  const lduMesh& mesh,
141  const labelList& value,
142  labelHashSet& sharedFaces
143 ) const
144 {
145  const lduAddressing& addr = mesh.lduAddr();
146  const labelUList& lower = addr.lowerAddr();
147  const labelUList& upper = addr.upperAddr();
148 
149  sharedFaces.clear();
150  sharedFaces.resize(addr.lowerAddr().size()/100);
151 
152  // Detect any faces inbetween same value
153  forAll(lower, facei)
154  {
155  label lowerData = value[lower[facei]];
156  label upperData = value[upper[facei]];
157 
158  if (lowerData != -1 && lowerData == upperData)
159  {
160  sharedFaces.insert(facei);
161  }
162  }
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
167 
168 Foam::MGridGenGAMGAgglomeration::MGridGenGAMGAgglomeration
169 (
170  const lduMesh& mesh,
171  const dictionary& controlDict
172 )
173 :
175  fvMesh_(refCast<const fvMesh>(mesh))
176 {
177  // Min, max size of agglomerated cells
178  label minSize(controlDict.get<label>("minSize"));
179  label maxSize(controlDict.get<label>("maxSize"));
180 
181  // Number of iterations applied to improve agglomeration consistency across
182  // processor boundaries
183  label nProcConsistencyIter
184  (
185  controlDict.get<label>("nProcConsistencyIter")
186  );
187 
188  // Start geometric agglomeration from the cell volumes and areas of the mesh
189  scalarField* VPtr = const_cast<scalarField*>(&fvMesh_.cellVolumes());
190 
191  scalarField magFaceAreas(sqrt(3.0)*mag(fvMesh_.faceAreas()));
192  SubField<scalar> magSf(magFaceAreas, fvMesh_.nInternalFaces());
193 
194  scalarField* magSfPtr = const_cast<scalarField*>
195  (
196  &magSf.operator const scalarField&()
197  );
198 
199  // Create the boundary area cell field
200  scalarField* magSbPtr(new scalarField(fvMesh_.nCells(), Zero));
201 
202  {
203  scalarField& magSb = *magSbPtr;
204 
205  const labelList& own = fvMesh_.faceOwner();
206  const vectorField& Sf = fvMesh_.faceAreas();
207 
208  forAll(Sf, facei)
209  {
210  if (!fvMesh_.isInternalFace(facei))
211  {
212  magSb[own[facei]] += mag(Sf[facei]);
213  }
214  }
215  }
216 
217  // Agglomerate until the required number of cells in the coarsest level
218  // is reached
219 
220  label nCreatedLevels = 0;
221 
222  while (nCreatedLevels < maxLevels_ - 1)
223  {
224  label nCoarseCells = -1;
225 
226  tmp<labelField> finalAgglomPtr = agglomerate
227  (
228  nCoarseCells,
229  minSize,
230  maxSize,
231  meshLevel(nCreatedLevels).lduAddr(),
232  *VPtr,
233  *magSfPtr,
234  *magSbPtr
235  );
236 
237  // Adjust weights only
238  for (int i=0; i<nProcConsistencyIter; i++)
239  {
240  const lduMesh& mesh = meshLevel(nCreatedLevels);
241  const lduAddressing& addr = mesh.lduAddr();
242  const lduInterfacePtrsList interfaces = mesh.interfaces();
243 
244  const labelField& agglom = finalAgglomPtr();
245 
246  // Global numbering
247  const globalIndex globalNumbering(nCoarseCells);
248 
249  labelField globalAgglom(addr.size());
250  forAll(agglom, celli)
251  {
252  globalAgglom[celli] = globalNumbering.toGlobal(agglom[celli]);
253  }
254 
255  // Get the interface agglomeration
256  PtrList<labelList> nbrGlobalAgglom;
257  swap(interfaces, globalAgglom, nbrGlobalAgglom);
258 
259 
260  // Get the interface agglomeration on a cell basis (-1 for all
261  // other cells)
262  labelList cellToNbrAgglom;
263  getNbrAgglom(addr, interfaces, nbrGlobalAgglom, cellToNbrAgglom);
264 
265 
266  // Mark all faces inbetween cells with same nbragglomeration
267  labelHashSet sharedFaces(addr.size()/100);
268  detectSharedFaces(mesh, cellToNbrAgglom, sharedFaces);
269 
270 
271  //- Note: in-place update of weights is more effective it seems?
272  // Should not be. fluke?
273  //scalarField weights(*faceWeightsPtr);
274  scalarField weights = *magSfPtr;
275  for (const label facei : sharedFaces)
276  {
277  weights[facei] *= 2.0;
278  }
279 
280  // Redo the agglomeration using the new weights
281  finalAgglomPtr = agglomerate
282  (
283  nCoarseCells,
284  minSize,
285  maxSize,
286  meshLevel(nCreatedLevels).lduAddr(),
287  *VPtr,
288  weights,
289  *magSbPtr
290  );
291  }
292 
293  if (continueAgglomerating(finalAgglomPtr().size(), nCoarseCells))
294  {
295  nCells_[nCreatedLevels] = nCoarseCells;
296  restrictAddressing_.set(nCreatedLevels, finalAgglomPtr);
297  }
298  else
299  {
300  break;
301  }
302 
303  agglomerateLduAddressing(nCreatedLevels);
304 
305  // Agglomerate the cell volumes field for the next level
306  {
307  scalarField* aggVPtr
308  (
309  new scalarField(meshLevels_[nCreatedLevels].size())
310  );
311 
312  // Restrict but no parallel agglomeration (not supported)
313  restrictField(*aggVPtr, *VPtr, nCreatedLevels, false);
314 
315  if (nCreatedLevels)
316  {
317  delete VPtr;
318  }
319 
320  VPtr = aggVPtr;
321  }
322 
323  // Agglomerate the face areas field for the next level
324  {
325  scalarField* aggMagSfPtr
326  (
327  new scalarField
328  (
329  meshLevels_[nCreatedLevels].upperAddr().size(),
330  0
331  )
332  );
333 
334  restrictFaceField(*aggMagSfPtr, *magSfPtr, nCreatedLevels);
335 
336  if (nCreatedLevels)
337  {
338  delete magSfPtr;
339  }
340 
341  magSfPtr = aggMagSfPtr;
342  }
343 
344  // Agglomerate the cell boundary areas field for the next level
345  {
346  scalarField* aggMagSbPtr
347  (
348  new scalarField(meshLevels_[nCreatedLevels].size())
349  );
350 
351  // Restrict but no parallel agglomeration (not supported)
352  restrictField(*aggMagSbPtr, *magSbPtr, nCreatedLevels, false);
353 
354  delete magSbPtr;
355  magSbPtr = aggMagSbPtr;
356  }
357 
358  nCreatedLevels++;
359  }
360 
361  // Shrink the storage of the levels to those created
362  compactLevels(nCreatedLevels);
363 
364  // Delete temporary geometry storage
365  if (nCreatedLevels)
366  {
367  delete VPtr;
368  delete magSfPtr;
369  }
370  delete magSbPtr;
371 }
372 
373 
374 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:114
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::GAMGAgglomeration
Geometric agglomerated algebraic multigrid agglomeration class.
Definition: GAMGAgglomeration.H:64
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::UPstream::waitRequests
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:262
Foam::HashSet< label, Hash< label > >
Foam::lduAddressing::size
label size() const
Return number of equations.
Definition: lduAddressing.H:171
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
processorLduInterface.H
Foam::SubField
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: Field.H:64
Foam::stringOps::lower
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Definition: stringOps.C:1184
Foam::Field< scalar >
Foam::fvMesh::interfaces
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:708
Foam::fvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:691
controlDict
runTime controlDict().readEntry("adjustTimeStep"
Definition: debug.C:143
Foam::UPtrList< const lduInterface >
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
MGridGenGAMGAgglomeration.H
Foam::UPstream::commsTypes::nonBlocking
Foam::UPstream::nRequests
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:252
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< label >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::stringOps::upper
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Definition: stringOps.C:1200
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::lduMesh
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:62
Foam::lduInterfacePtrsList
UPtrList< const lduInterface > lduInterfacePtrsList
List of coupled interface fields to be used in coupling.
Definition: lduInterfacePtrsList.H:44
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:240