pairPatchAgglomeration.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-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "meshTools.H"
31#include "edgeHashes.H"
32#include "unitConversion.H"
33
34// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35
36void Foam::pairPatchAgglomeration::compactLevels(const label nCreatedLevels)
37{
38 nFaces_.setSize(nCreatedLevels);
39 restrictAddressing_.setSize(nCreatedLevels);
40 patchLevels_.setSize(nCreatedLevels);
41}
42
43
44bool Foam::pairPatchAgglomeration::continueAgglomerating
45(
46 const label nLocal,
47 const label nLocalOld
48)
49{
50 // Keep agglomerating
51 // - if global number of faces is still changing
52 // - and if local number of faces still too large (on any processor)
53 // or if global number of faces still too large
54
55 label nGlobal = returnReduce(nLocal, sumOp<label>());
56 label nGlobalOld = returnReduce(nLocalOld, sumOp<label>());
57
58 return
59 (
60 returnReduce(nLocal > nFacesInCoarsestLevel_, orOp<bool>())
61 || nGlobal > nGlobalFacesInCoarsestLevel_
62 )
63 && nGlobal != nGlobalOld;
64}
65
66
67void Foam::pairPatchAgglomeration::setLevel0EdgeWeights()
68{
69 const bPatch& coarsePatch = patchLevels_[0];
70 forAll(coarsePatch.edges(), i)
71 {
72 if (coarsePatch.isInternalEdge(i))
73 {
74 scalar edgeLength =
75 coarsePatch.edges()[i].mag(coarsePatch.localPoints());
76
77 const labelList& eFaces = coarsePatch.edgeFaces()[i];
78
79 if (eFaces.size() == 2)
80 {
81 scalar cosI =
82 coarsePatch.faceNormals()[eFaces[0]]
83 & coarsePatch.faceNormals()[eFaces[1]];
84
85 const edge edgeCommon = edge(eFaces[0], eFaces[1]);
86
87 if (facePairWeight_.found(edgeCommon))
88 {
89 facePairWeight_[edgeCommon] += edgeLength;
90 }
91 else
92 {
93 facePairWeight_.insert(edgeCommon, edgeLength);
94 }
95
96 if (mag(cosI) < Foam::cos(degToRad(featureAngle_)))
97 {
98 facePairWeight_[edgeCommon] = -1.0;
99 }
100 }
101 else
102 {
103 forAll(eFaces, j)
104 {
105 for (label k = j+1; k<eFaces.size(); k++)
106 {
107 facePairWeight_.insert
108 (
109 edge(eFaces[j], eFaces[k]),
110 -1.0
111 );
112 }
113 }
114 }
115 }
116 }
117}
118
119
120void Foam::pairPatchAgglomeration::setEdgeWeights
121(
122 const label fineLevelIndex
123)
124{
125 const bPatch& coarsePatch = patchLevels_[fineLevelIndex];
126 const labelList& fineToCoarse = restrictAddressing_[fineLevelIndex];
127 const label nCoarseI = max(fineToCoarse) + 1;
128 labelListList coarseToFine(invertOneToMany(nCoarseI, fineToCoarse));
129
130 edgeHashSet fineFeaturedFaces(coarsePatch.nEdges()/10);
131
132 // Map fine faces with featured edge into coarse faces
133 forAllConstIters(facePairWeight_, iter)
134 {
135 if (iter() == -1.0)
136 {
137 const edge e = iter.key();
138 const edge edgeFeatured
139 (
140 fineToCoarse[e[0]],
141 fineToCoarse[e[1]]
142 );
143 fineFeaturedFaces.insert(edgeFeatured);
144 }
145 }
146
147 // Clean old weights
148 facePairWeight_.clear();
149 facePairWeight_.resize(coarsePatch.nEdges());
150
151 forAll(coarsePatch.edges(), i)
152 {
153 if (coarsePatch.isInternalEdge(i))
154 {
155 scalar edgeLength =
156 coarsePatch.edges()[i].mag(coarsePatch.localPoints());
157
158 const labelList& eFaces = coarsePatch.edgeFaces()[i];
159
160 if (eFaces.size() == 2)
161 {
162 const edge edgeCommon = edge(eFaces[0], eFaces[1]);
163 if (facePairWeight_.found(edgeCommon))
164 {
165 facePairWeight_[edgeCommon] += edgeLength;
166 }
167 else
168 {
169 facePairWeight_.insert(edgeCommon, edgeLength);
170 }
171 // If the fine 'pair' faces was featured edge so it is
172 // the coarse 'pair'
173 if (fineFeaturedFaces.found(edgeCommon))
174 {
175 facePairWeight_[edgeCommon] = -1.0;
176 }
177 }
178 else
179 {
180 // Set edge as barrier by setting weight to -1
181 forAll(eFaces, j)
182 {
183 for (label k = j+1; k<eFaces.size(); k++)
184 {
185 facePairWeight_.insert
186 (
187 edge(eFaces[j], eFaces[k]),
188 -1.0
189 );
190 }
191 }
192 }
193 }
194 }
195}
196
197
198// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
199
201(
202 const faceList& faces,
203 const pointField& points,
205)
206:
207 mergeLevels_
208 (
209 controlDict.getOrDefault<label>("mergeLevels", 2)
210 ),
211 maxLevels_(50),
212 nFacesInCoarsestLevel_
213 (
214 controlDict.get<label>("nFacesInCoarsestLevel")
215 ),
216 nGlobalFacesInCoarsestLevel_(labelMax),
217 //(
218 // controlDict.get<label>("nGlobalFacesInCoarsestLevel")
219 //),
220 featureAngle_
221 (
222 controlDict.getOrDefault<scalar>("featureAngle", 0)
223 ),
224 nFaces_(maxLevels_),
225 restrictAddressing_(maxLevels_),
226 restrictTopBottomAddressing_(identity(faces.size())),
227 patchLevels_(maxLevels_),
228 facePairWeight_(faces.size())
229{
230 // Set base fine patch
231 patchLevels_.set(0, new bPatch(faces, points));
232
233 // Set number of faces for the base patch
234 nFaces_[0] = faces.size();
235
236 // Set edge weights for level 0
237 setLevel0EdgeWeights();
238}
239
240
242(
243 const faceList& faces,
244 const pointField& points,
245 const label mergeLevels,
246 const label maxLevels,
247 const label nFacesInCoarsestLevel, // local number of cells
248 const label nGlobalFacesInCoarsestLevel, // global number of cells
249 const scalar featureAngle
250)
251:
252 mergeLevels_(mergeLevels),
253 maxLevels_(maxLevels),
254 nFacesInCoarsestLevel_(nFacesInCoarsestLevel),
255 nGlobalFacesInCoarsestLevel_(nGlobalFacesInCoarsestLevel),
256 featureAngle_(featureAngle),
257 nFaces_(maxLevels_),
258 restrictAddressing_(maxLevels_),
259 restrictTopBottomAddressing_(identity(faces.size())),
260 patchLevels_(maxLevels_),
261 facePairWeight_(faces.size())
262{
263 // Set base fine patch
264 patchLevels_.set(0, new bPatch(faces, points));
265
266 // Set number of faces for the base patch
267 nFaces_[0] = faces.size();
268
269 // Set edge weights for level 0
270 setLevel0EdgeWeights();
271}
272
273
274// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
275
277{}
278
279
280// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
281
284(
285 const label i
286) const
287{
288 return patchLevels_[i];
289}
290
291
292void Foam::pairPatchAgglomeration::mapBaseToTopAgglom
293(
294 const label fineLevelIndex
295)
296{
297 const labelList& fineToCoarse = restrictAddressing_[fineLevelIndex];
298 forAll(restrictTopBottomAddressing_, i)
299 {
300 restrictTopBottomAddressing_[i] =
301 fineToCoarse[restrictTopBottomAddressing_[i]];
302 }
303}
304
305
306bool Foam::pairPatchAgglomeration::agglomeratePatch
307(
308 const bPatch& patch,
309 const labelList& fineToCoarse,
310 const label fineLevelIndex
311)
312{
313 if (min(fineToCoarse) == -1)
314 {
316 << "min(fineToCoarse) == -1" << exit(FatalError);
317 }
318
319 if (fineToCoarse.size() == 0)
320 {
321 return false;
322 }
323
324 if (fineToCoarse.size() != patch.size())
325 {
327 << "restrict map does not correspond to fine level. " << endl
328 << " Sizes: restrictMap: " << fineToCoarse.size()
329 << " nEqns: " << patch.size()
330 << abort(FatalError);
331 }
332
333 const label nCoarseI = max(fineToCoarse) + 1;
334 List<face> patchFaces(nCoarseI);
335
336
337 // Patch faces per agglomeration
338 labelListList coarseToFine(invertOneToMany(nCoarseI, fineToCoarse));
339
340 for (label coarseI = 0; coarseI < nCoarseI; coarseI++)
341 {
342 const labelList& fineFaces = coarseToFine[coarseI];
343
344 // Construct single face
346 (
347 IndirectList<face>(patch, fineFaces),
348 patch.points()
349 );
350
351 if (upp.edgeLoops().size() != 1)
352 {
353 if (fineFaces.size() == 2)
354 {
355 const edge e(fineFaces[0], fineFaces[1]);
356 facePairWeight_[e] = -1.0;
357 }
358 else if (fineFaces.size() == 3)
359 {
360 const edge e(fineFaces[0], fineFaces[1]);
361 const edge e1(fineFaces[0], fineFaces[2]);
362 const edge e2(fineFaces[2], fineFaces[1]);
363 facePairWeight_[e] = -1.0;
364 facePairWeight_[e1] = -1.0;
365 facePairWeight_[e2] = -1.0;
366 }
367 return false;
368 }
369
370 patchFaces[coarseI] = face
371 (
373 (
374 upp.meshPoints(),
375 upp.edgeLoops()[0]
376 )
377 );
378 }
379
380 patchLevels_.set
381 (
382 fineLevelIndex,
383 new bPatch
384 (
385 SubList<face>(patchFaces, nCoarseI, 0),
386 patch.points()
387 )
388 );
389
390 return true;
391}
392
393
395{
396 label nPairLevels = 0;
397 label nCreatedLevels = 1; // 0 level is the base patch
398
399 label nCoarseFaces = 0;
400 label nCoarseFacesOld = 0;
401
402 while (nCreatedLevels < maxLevels_)
403 {
404 const bPatch& patch = patchLevels_[nCreatedLevels - 1];
405
406 // Agglomerate locally
407 tmp<labelField> tfinalAgglom;
408
409 bool createdLevel = false;
410 while (!createdLevel)
411 {
412 // Agglomerate locally using edge weights
413 // - calculates nCoarseFaces; returns fine to coarse addressing
414 tfinalAgglom = agglomerateOneLevel(nCoarseFaces, patch);
415
416 if (nCoarseFaces == 0)
417 {
418 break;
419 }
420 else
421 {
422 // Attempt to create coarse face addressing
423 // - returns true if successful; otherwise resets edge weights
424 // and assume try again...
425 createdLevel = agglomeratePatch
426 (
427 patch,
428 tfinalAgglom,
429 nCreatedLevels
430 );
431 }
432 }
433
434 if (createdLevel)
435 {
436 restrictAddressing_.set(nCreatedLevels, tfinalAgglom);
437
438 mapBaseToTopAgglom(nCreatedLevels);
439
440 setEdgeWeights(nCreatedLevels);
441
442 if (nPairLevels % mergeLevels_)
443 {
444 combineLevels(nCreatedLevels);
445 }
446 else
447 {
448 nCreatedLevels++;
449 }
450
451 nPairLevels++;
452
453 nFaces_[nCreatedLevels] = nCoarseFaces;
454 }
455
456 // Check to see if we need to continue agglomerating
457 // - Note: performs parallel reductions
458 if (!continueAgglomerating(nCoarseFaces, nCoarseFacesOld))
459 {
460 break;
461 }
462
463 nCoarseFacesOld = nCoarseFaces;
464 }
465
466 compactLevels(nCreatedLevels);
467}
468
469
470Foam::tmp<Foam::labelField> Foam::pairPatchAgglomeration::agglomerateOneLevel
471(
472 label& nCoarseFaces,
473 const bPatch& patch
474)
475{
476 const label nFineFaces = patch.size();
477
478 tmp<labelField> tcoarseCellMap(new labelField(nFineFaces, -1));
479 labelField& coarseCellMap = tcoarseCellMap.ref();
480
481 const labelListList& faceFaces = patch.faceFaces();
482
483 nCoarseFaces = 0;
484
485 forAll(faceFaces, facei)
486 {
487 const labelList& fFaces = faceFaces[facei];
488
489 if (coarseCellMap[facei] < 0)
490 {
491 label matchFaceNo = -1;
492 label matchFaceNeibNo = -1;
493 scalar maxFaceWeight = -GREAT;
494
495 // Check faces to find ungrouped neighbour with largest face weight
496 forAll(fFaces, i)
497 {
498 label faceNeig = fFaces[i];
499 const edge edgeCommon = edge(facei, faceNeig);
500 if
501 (
502 facePairWeight_[edgeCommon] > maxFaceWeight
503 && coarseCellMap[faceNeig] < 0
504 && facePairWeight_[edgeCommon] != -1.0
505 )
506 {
507 // Match found. Pick up all the necessary data
508 matchFaceNo = facei;
509 matchFaceNeibNo = faceNeig;
510 maxFaceWeight = facePairWeight_[edgeCommon];
511 }
512 }
513
514 if (matchFaceNo >= 0)
515 {
516 // Make a new group
517 coarseCellMap[matchFaceNo] = nCoarseFaces;
518 coarseCellMap[matchFaceNeibNo] = nCoarseFaces;
519 nCoarseFaces++;
520 }
521 else
522 {
523 // No match. Find the best neighbouring cluster and
524 // put the cell there
525 label clusterMatchFaceNo = -1;
526 scalar clusterMaxFaceCoeff = -GREAT;
527
528 forAll(fFaces, i)
529 {
530 label faceNeig = fFaces[i];
531 const edge edgeCommon = edge(facei, faceNeig);
532 if
533 (
534 facePairWeight_[edgeCommon] > clusterMaxFaceCoeff
535 && facePairWeight_[edgeCommon] != -1.0
536 && coarseCellMap[faceNeig] >= 0
537 )
538 {
539 clusterMatchFaceNo = faceNeig;
540 clusterMaxFaceCoeff = facePairWeight_[edgeCommon];
541 }
542 }
543
544 if (clusterMatchFaceNo > 0)
545 {
546 // Add the cell to the best cluster
547 coarseCellMap[facei] = coarseCellMap[clusterMatchFaceNo];
548 }
549 else
550 {
551 // If not create single-cell "clusters" for each
552 coarseCellMap[facei] = nCoarseFaces;
553 nCoarseFaces++;
554 }
555 }
556 }
557 }
558
559 // Check that all faces are part of clusters,
560 for (label facei=0; facei<nFineFaces; facei++)
561 {
562 if (coarseCellMap[facei] < 0)
563 {
565 << " face " << facei
566 << " is not part of a cluster"
567 << exit(FatalError);
568 }
569 }
570
571 return tcoarseCellMap;
572}
573
574
575void Foam::pairPatchAgglomeration::combineLevels(const label curLevel)
576{
577 label prevLevel = curLevel - 1;
578
579 // Set the previous level nCells to the current
580 nFaces_[prevLevel] = nFaces_[curLevel];
581
582 // Map the restrictAddressing from the coarser level into the previous
583 // finer level
584
585 const labelList& curResAddr = restrictAddressing_[curLevel];
586 labelList& prevResAddr = restrictAddressing_[prevLevel];
587
588 forAll(prevResAddr, i)
589 {
590 prevResAddr[i] = curResAddr[prevResAddr[i]];
591 }
592
593 // Delete the restrictAddressing for the coarser level
594 restrictAddressing_.set(curLevel, nullptr);
595
596 patchLevels_.set(prevLevel, patchLevels_.set(curLevel, nullptr));
597}
598
599
600// ************************************************************************* //
label k
void setSize(const label n)
Alias for resize()
Definition: List.H:218
A list of faces which address into the list of points.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Primitive patch pair agglomerate method.
PtrList< bPatch > patchLevels_
Hierarchy of patch addressing.
labelList nFaces_
The number of faces in each level.
const bPatch & patchLevel(const label leveli) const
Return primitivePatch of given level.
void agglomerate()
Agglomerate patch.
PtrList< labelField > restrictAddressing_
Cell restriction addressing array.
PrimitivePatch< List< face >, const pointField > bPatch
A class for managing temporary objects.
Definition: tmp.H:65
runTime controlDict().readEntry("adjustTimeStep"
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
const std::string patch
OpenFOAM patch number as a std::string.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
List< label > labelList
A List of labels.
Definition: List.H:66
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values (not the indices) of a list.
constexpr label labelMax
Definition: label.H:61
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
HashSet< edge, Hash< edge > > edgeHashSet
A HashSet with edge for its key.
Definition: edgeHashes.H:51
PrimitivePatch<::Foam::List< face >, const pointField > bPatch
error FatalError
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:114
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensionedScalar cos(const dimensionedScalar &ds)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Unit conversion functions.