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