holeToFace.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) 2020-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 "holeToFace.H"
29#include "transform.H"
30#include "faceSet.H"
31#include "cellSet.H"
33#include "OBJstream.H"
34//#include "fvMesh.H"
35//#include "volFields.H"
36//#include "surfaceFields.H"
37#include "topoDistanceData.H"
38#include "FaceCellWave.H"
39#include "syncTools.H"
40
42#include "PatchEdgeFaceWave.H"
44
45// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46
47namespace Foam
48{
55 (
58 word,
59 hole
60 );
62 (
65 istream,
66 hole
67 );
68}
69
70
71Foam::topoSetSource::addToUsageTable Foam::holeToFace::usage_
72(
73 holeToFace::typeName,
74 "\n Usage: holeToFace <faceSet> ((x0 y0 z0) (x1 y1 z1))\n\n"
75 " Select faces disconnecting the individual regions"
76 " (each indicated by a locations).\n"
77);
78
79
80// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
81
82//Foam::label Foam::holeToFace::globalCount
83//(
84// const bitSet& isMasterFace,
85// const bitSet& set
86//)
87//{
88// if (set.size() != isMasterFace.size())
89// {
90// FatalErrorInFunction << "problem" << exit(FatalError);
91// }
92//
93// label n = 0;
94// for (const label facei : set)
95// {
96// if (isMasterFace(facei))
97// {
98// n++;
99// }
100// }
101// return returnReduce(n, sumOp<label>());
102//}
103
104
105//void Foam::holeToFace::checkFaceSync
106//(
107// const string& setName,
108// const bitSet& set
109//) const
110//{
111// if (set.size() != mesh_.nFaces())
112// {
113// FatalErrorInFunction<< "problem" << exit(FatalError);
114// }
115// bitSet orSet(set);
116// syncTools::syncFaceList(mesh_, orSet, orEqOp<unsigned int>());
117// bitSet andSet(set);
118// syncTools::syncFaceList(mesh_, andSet, andEqOp<unsigned int>());
119//
120// forAll(orSet, facei)
121// {
122// if (orSet[facei] != andSet[facei] || orSet[facei] != set[facei])
123// {
124// FatalErrorInFunction<< "problem for set " << setName
125// << "face:" << facei
126// << " at:" << mesh_.faceCentres()[facei]
127// << " patch:" << mesh_.boundaryMesh().whichPatch(facei)
128// << " set:" << set[facei]
129// << " orSet:" << orSet[facei]
130// << " andSet:" << andSet[facei]
131// << exit(FatalError);
132// }
133// }
134//}
135
136
137//void Foam::holeToFace::checkFaceSync
138//(
139// const string& fldName,
140// const labelList& fld
141//) const
142//{
143// if (fld.size() != mesh_.nFaces())
144// {
145// FatalErrorInFunction<< "problem" << exit(FatalError);
146// }
147// labelList maxFld(fld);
148// syncTools::syncFaceList(mesh_, maxFld, maxEqOp<label>());
149// labelList minFld(fld);
150// syncTools::syncFaceList(mesh_, minFld, minEqOp<label>());
151//
152// forAll(maxFld, facei)
153// {
154// if (maxFld[facei] != minFld[facei] || maxFld[facei] != fld[facei])
155// {
156// FatalErrorInFunction<< "problem for field " << fldName
157// << "face:" << facei
158// << " at:" << mesh_.faceCentres()[facei]
159// << " patch:" << mesh_.boundaryMesh().whichPatch(facei)
160// << " fld:" << fld[facei]
161// << " maxFld:" << maxFld[facei]
162// << " minFld:" << minFld[facei]
163// << exit(FatalError);
164// }
165// }
166//}
167
168
169//void Foam::holeToFace::checkFaceSync
170//(
171// const string& fldName,
172// const List<unsigned int>& fld
173//) const
174//{
175// if (fld.size() != mesh_.nFaces())
176// {
177// FatalErrorInFunction<< "problem" << exit(FatalError);
178// }
179// List<unsigned int> orFld(fld);
180// syncTools::syncFaceList(mesh_, orFld, bitOrEqOp<unsigned int>());
181// List<unsigned int> andFld(fld);
182// syncTools::syncFaceList(mesh_, andFld, bitAndEqOp<unsigned int>());
183// forAll(orFld, facei)
184// {
185// if (orFld[facei] != andFld[facei] || orFld[facei] != fld[facei])
186// {
187// FatalErrorInFunction<< "problem for field " << fldName
188// << "face:" << facei
189// << " at:" << mesh_.faceCentres()[facei]
190// << " patch:" << mesh_.boundaryMesh().whichPatch(facei)
191// << " fld:" << fld[facei]
192// << " orFld:" << orFld[facei]
193// << " andFld:" << andFld[facei]
194// << exit(FatalError);
195// }
196// }
197//}
198
199
200//void Foam::holeToFace::writeCellField
201//(
202// const word& name,
203// const labelList& labelFld
204//) const
205//{
206// Pout<< "Writing field " << name << endl;
207// if (labelFld.size() != mesh_.nCells())
208// {
209// FatalErrorInFunction << exit(FatalError);
210// }
211//
212// const fvMesh& fvm = dynamic_cast<const fvMesh&>(mesh_);
213//
214// volScalarField fld
215// (
216// IOobject
217// (
218// name,
219// fvm.time().timeName(),
220// fvm,
221// IOobject::NO_READ,
222// IOobject::NO_WRITE,
223// false
224// ),
225// fvm,
226// dimensionedScalar("zero", dimless, scalar(0))
227// );
228// forAll(labelFld, i)
229// {
230// fld[i] = labelFld[i];
231// }
232// fld.correctBoundaryConditions();
233// fld.write();
234//}
235
236
237//void Foam::holeToFace::writeFaceField
238//(
239// const word& name,
240// const labelList& labelFld
241//) const
242//{
243// Pout<< "Writing field " << name << endl;
244// if (labelFld.size() != mesh_.nFaces())
245// {
246// FatalErrorInFunction << exit(FatalError);
247// }
248//
249// const fvMesh& fvm = dynamic_cast<const fvMesh&>(mesh_);
250//
251// surfaceScalarField fld
252// (
253// IOobject
254// (
255// name,
256// fvm.time().timeName(),
257// fvm,
258// IOobject::NO_READ,
259// IOobject::NO_WRITE,
260// false
261// ),
262// fvm,
263// dimensionedScalar("zero", dimless, scalar(0))
264// );
265// for (label i = 0; i < mesh_.nInternalFaces(); i++)
266// {
267// fld[i] = labelFld[i];
268// }
269// surfaceScalarField::Boundary& bfld = fld.boundaryFieldRef();
270// forAll(bfld, patchi)
271// {
272// fvsPatchScalarField& pfld = bfld[patchi];
273// forAll(pfld, i)
274// {
275// pfld[i] = labelFld[pfld.patch().start()+i];
276// }
277// }
278// fld.write();
279//
280// // Write as faceSet as well
281// const bitSet isMasterFace(syncTools::getInternalOrMasterFaces(mesh_));
282//
283// faceSet set(mesh_, name, 100);
284// label nMasters = 0;
285// forAll(labelFld, facei)
286// {
287// if (labelFld[facei] >= 0)
288// {
289// set.insert(facei);
290// if (isMasterFace(facei))
291// {
292// nMasters++;
293// }
294// }
295// }
296// Pout<< "Writing " << returnReduce(nMasters, sumOp<label>())
297// << " >= 0 faces to faceSet " << set.name() << endl;
298// set.write();
299//}
300
301
302void Foam::holeToFace::writeFaces
303(
304 const word& name,
305 const bitSet& faceFld
306) const
307{
309 OBJstream str(mesh_.time().timePath()/name);
310 Pout<< "Writing " << faceFld.count() << " faces to " << str.name() << endl;
311
312 for (const label facei : faceFld)
313 {
314 str.write(mesh_.faces()[facei], mesh_.points(), false);
315 }
316}
317
318
319void Foam::holeToFace::calculateDistance
320(
321 const labelList& seedFaces,
322 const bitSet& isBlockedCell,
323 const bitSet& isBlockedFace,
324 labelList& cellDist,
325 labelList& faceDist
326) const
327{
328 if (isBlockedCell.size() != mesh_.nCells())
329 {
330 FatalErrorInFunction << "Problem" << exit(FatalError);
331 }
332 if (isBlockedFace.size() != mesh_.nFaces())
333 {
334 FatalErrorInFunction << "Problem" << exit(FatalError);
335 }
336
337 //const bitSet isMasterFace(syncTools::getInternalOrMasterFaces(mesh_));
338
339 // Field on cells and faces.
340 List<topoDistanceData<label>> cellData(mesh_.nCells());
341 List<topoDistanceData<label>> faceData(mesh_.nFaces());
342
343 // Start of changes
344 List<topoDistanceData<label>> seedData
345 (
346 seedFaces.size(),
347 topoDistanceData<label>(0, 123)
348 );
349 //Pout<< "Seeded "
350 // << globalCount(isMasterFace, bitSet(mesh_.nFaces(), seedFaces))
351 // << " out of " << returnReduce(mesh_.nFaces(), sumOp<label>()) << endl;
352
353 // Make sure we don't walk through inactive cells
354 //Pout<< "blocking "
355 // << returnReduce(isBlockedCell.count(), sumOp<label>())
356 // << " cells" << endl;
357 for (const label celli : isBlockedCell)
358 {
359 cellData[celli] = topoDistanceData<label>(0, 0);
360 }
361 //Pout<< "blocking "
362 // << globalCount(isMasterFace, isBlockedFace) << " faces" << endl;
363 for (const label facei : isBlockedFace)
364 {
365 faceData[facei] = topoDistanceData<label>(0, 0);
366 }
367
368 // Propagate information inwards
369 FaceCellWave<topoDistanceData<label>> deltaCalc
370 (
371 mesh_,
372 seedFaces,
373 seedData,
374 faceData,
375 cellData,
376 mesh_.globalData().nTotalCells()+1
377 );
378
379 // And extract
380 //bool haveWarned = false;
381 forAll(cellData, celli)
382 {
383 if (!isBlockedCell[celli])
384 {
385 if (!cellData[celli].valid(deltaCalc.data()))
386 {
387 //if (!haveWarned)
388 //{
389 // WarningInFunction
390 // << "Did not visit some cells, e.g. cell " << celli
391 // << " at " << mesh_.cellCentres()[celli] << endl;
392 // haveWarned = true;
393 //}
394 }
395 else
396 {
397 cellDist[celli] = cellData[celli].distance();
398 }
399 }
400 }
401
402 forAll(faceDist, facei)
403 {
404 if (!isBlockedFace[facei])
405 {
406 if (!faceData[facei].valid(deltaCalc.data()))
407 {
408 //if (!haveWarned)
409 //{
410 // WarningInFunction
411 // << "Did not visit some faces, e.g. face " << facei
412 // << " at " << mesh_.faceCentres()[facei] << endl;
413 // haveWarned = true;
414 //}
415 }
416 else
417 {
418 faceDist[facei] = faceData[facei].distance();
419 }
420 }
421 }
422}
423
424
425Foam::bitSet Foam::holeToFace::frontFaces
426(
427 const bitSet& isSurfaceFace,
428 const List<unsigned int>& locationFaces,
429 const bitSet& isHoleCell
430) const
431{
432 const labelList& faceOwner = mesh_.faceOwner();
433 const labelList& faceNeighbour = mesh_.faceNeighbour();
434 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
435
436 bitSet isFrontFace(mesh_.nFaces());
437 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
438 {
439 if (!isSurfaceFace[facei])
440 {
441 const label ownHole = isHoleCell[faceOwner[facei]];
442 const label neiHole = isHoleCell[faceNeighbour[facei]];
443
444 if (ownHole != neiHole)
445 {
446 unsigned int masks = locationFaces[facei];
447 if (masks == 0u)
448 {
449 FatalErrorInFunction << "face:" << facei
450 << " at:" << mesh_.faceCentres()[facei]
451 << " not any front" << exit(FatalError);
452 }
453
454 // Count number of bits set
455 const label nSet = BitOps::bit_count(masks);
456
457 if (nSet == 1)
458 {
459 isFrontFace.set(facei);
460 }
461 }
462 }
463 }
464
465 // Get neighbouring cell data
466 bitSet isHoleNeiCell(mesh_.nBoundaryFaces());
467 {
468 for (const polyPatch& pp : pbm)
469 {
470 label bFacei = pp.start()-mesh_.nInternalFaces();
471 const labelUList& faceCells = pp.faceCells();
472
473 for (const label celli : faceCells)
474 {
475 isHoleNeiCell[bFacei] = isHoleCell[celli];
476 ++bFacei;
477 }
478 }
479 syncTools::swapBoundaryFaceList(mesh_, isHoleNeiCell);
480 }
481
482
483 for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
484 {
485 if (!isSurfaceFace[facei])
486 {
487 const label ownHole = isHoleCell[faceOwner[facei]];
488 const label neiHole = isHoleNeiCell[facei-mesh_.nInternalFaces()];
489
490 if (ownHole != neiHole)
491 {
492 unsigned int masks = locationFaces[facei];
493 if (masks == 0u)
494 {
495 FatalErrorInFunction << "face:" << facei
496 << " at:" << mesh_.faceCentres()[facei]
497 << " not any front" << exit(FatalError);
498 }
499
500 // Count number of bits set
501 const label nSet = BitOps::bit_count(masks);
502
503 if (nSet == 1)
504 {
505 isFrontFace.set(facei);
506 }
507 }
508 }
509 }
510 return isFrontFace;
511}
512
513
514Foam::bitSet Foam::holeToFace::findClosure
515(
516 const bitSet& isSurfaceFace, // intersected faces
517 const bitSet& isCandidateHoleCell, // starting blockage
518 const labelListList& locationCells // cells per zone
519) const
520{
521 const labelList& faceOwner = mesh_.faceOwner();
522 const labelList& faceNeighbour = mesh_.faceNeighbour();
523
524 if (zonePoints_.size() < 2)
525 {
526 FatalErrorInFunction << "single region only : "
527 << flatOutput(zonePoints_) << exit(FatalError);
528 }
529
530 if (zonePoints_.size() > 31)
531 {
532 FatalErrorInFunction << "only support < 32 locations in mesh."
533 << " Currently : " << flatOutput(zonePoints_) << exit(FatalError);
534 }
535
536
537 //const bitSet isMasterFace(syncTools::getInternalOrMasterFaces(mesh_));
538
539 bitSet isHoleCell(isCandidateHoleCell);
540 for (const labelList& zoneCells : locationCells)
541 {
542 for (const label celli : zoneCells)
543 {
544 if (celli != -1)
545 {
546 isHoleCell.unset(celli);
547 }
548 }
549 }
550
551 bitSet notHoleCell(isHoleCell);
552 notHoleCell.flip();
553
554 if (debug)
555 {
556 Pout<< "holeToFace::findClosure :"
557 << " locationCells:" << flatOutput(locationCells) << nl
558 << "holeToFace::findClosure :"
559 << " initial blocked faces:" << isSurfaceFace.count()
560 << " candidate closure cells:" << isHoleCell.count()
561 << endl;
562 }
563
564 // Distance to surface for every cell/face inside isHoleCell
565 labelList surfaceCellDist(mesh_.nCells(), -1);
566 labelList surfaceNeiCellDist(mesh_.nBoundaryFaces(), -1);
567 labelList surfaceFaceDist(mesh_.nFaces(), -1);
568 {
569 calculateDistance
570 (
571 isSurfaceFace.toc(), // seed faces
572 notHoleCell, // no need to walk through non-hole cells
573 bitSet(mesh_.nFaces()), // no blocked faces
574 surfaceCellDist,
575 surfaceFaceDist
576 );
578 (
579 mesh_,
580 surfaceCellDist,
581 surfaceNeiCellDist
582 );
583 //writeCellField("surfaceCellDistance0", surfaceCellDist);
584 //writeFaceField("surfaceFaceDistance0", surfaceFaceDist);
585 //checkFaceSync("surfaceFaceDistance0", surfaceFaceDist);
586 }
587
588 if (debug)
589 {
590 Pout<< "holeToFace::findClosure :"
591 << " calculated topological distance to initial blocked faces."
592 << " max distance:" << gMax(surfaceCellDist)
593 << endl;
594 }
595
596
597 // Find faces reachable from locationCells. If the locationCell is inside
598 // the blockage it will be only the faces of the cell. If it is outside
599 // it does a full walk to find the reachable faces on the outside of
600 // the blockage.
601 // Note: actual distance itself does not matter - only if they have been
602 // visited.
603 List<unsigned int> locationFaces(mesh_.nFaces(), 0u);
604 forAll(locationCells, zonei)
605 {
606 labelList cellDist(mesh_.nCells(), -1);
607 labelList faceDist(mesh_.nFaces(), -1);
608
609 labelList seedFaces;
610
611 const labelList& zoneCells = locationCells[zonei];
612 for (const label celli : zoneCells)
613 {
614 if (celli != -1)
615 {
616 cellDist[celli] = 0;
617 const cell& cFaces = mesh_.cells()[celli];
618 seedFaces.append(cFaces);
619 UIndirectList<label>(faceDist, cFaces) = 0;
620 }
621 }
622
623 // Extra par sync. Is this needed?
624 {
625 bitSet isSeedFace(mesh_.nFaces(), seedFaces);
627 (
628 mesh_,
629 isSeedFace,
630 orEqOp<unsigned int>()
631 );
632 seedFaces = isSeedFace.toc();
633 }
634
635
636 calculateDistance
637 (
638 seedFaces, // seed faces
639 isHoleCell, // do not walk through blocking cells
640 isSurfaceFace, // do not walk through surface
641 cellDist,
642 faceDist
643 );
644 //writeCellField("cellDistance" + Foam::name(zonei), cellDist);
645 //writeFaceField("faceDistance" + Foam::name(zonei), faceDist);
646 //checkFaceSync("faceDistance", faceDist);
647
648 // Add all reached faces
649 const unsigned int mask = (1u << zonei);
650 forAll(faceDist, facei)
651 {
652 if (faceDist[facei] >= 0)
653 {
654 locationFaces[facei] |= mask;
655 }
656 }
657 }
658
659 if (debug)
660 {
661 writeFaces("isSurfaceFace.obj", isSurfaceFace);
662 }
663 //if (debug)
664 //{
665 // bitSet isMultiBitFace(mesh_.nFaces());
666 // forAll(locationFaces, facei)
667 // {
668 // const unsigned int bits = locationFaces[facei];
669 // const label nSet = BitOps::bit_count(bits);
670 // if (nSet > 1)
671 // {
672 // isMultiBitFace.set(facei);
673 // }
674 // }
675 // writeFaces("isMultiBitFace.obj", isMultiBitFace);
676 //}
677
678
680 //- NOT CORRECT!!! At this point the walking has only done the distance
681 // outside the initial set of blocked faces. We'd have to
682 // walk through all faces before we can determine.
683 //bool haveLeak = false;
684 //forAll(locationFaces, facei)
685 //{
686 // if (!isSurfaceFace[facei])
687 // {
688 // const unsigned int bits = locationFaces[facei];
689 // const label nSet = BitOps::bit_count(bits);
690 // if (nSet > 1)
691 // {
692 // haveLeak = true;
693 //
694 // if (debug)
695 // {
696 // // Collect points
697 // DynamicList<pointField> connected;
698 // forAll(zonePoints_, zonei)
699 // {
700 // if (bits & (1u << zonei))
701 // {
702 // connected.append(zonePoints_[zonei]);
703 // }
704 // }
705 // Pout<< "holeToFace::findClosure :"
706 // << " found initial leak at face "
707 // << mesh_.faceCentres()[facei]
708 // << " between zones " << flatOutput(connected)
709 // << endl;
710 // }
711 // break;
712 // }
713 // }
714 //}
715 //reduce(haveLeak, orOp<bool>());
716 //
717 //if (!haveLeak)
718 //{
719 // if (debug)
720 // {
721 // Pout<< "holeToFace::findClosure :"
722 // << " did not find leak between zones "
723 // << flatOutput(zonePoints_) << endl;
724 // }
725 // return bitSet(mesh_.nFaces());
726 //}
727
728
729 // Front (on outside of hole cell but not connecting multiple locations
730 bitSet isFrontFace(frontFaces(isSurfaceFace, locationFaces, isHoleCell));
731
732
733 // Start off eroding the cells furthest away from the surface
734 label surfaceDist = gMax(surfaceCellDist);
735
736 // Work storage
737 //List<unsigned int> newLocationFaces(mesh_.nFaces());
738 //bitSet newHoleCell(mesh_.nCells());
739
740 while (surfaceDist >= 0)
741 {
742 // Erode cells with >= surfaceDist:
743 // - unmark cell as blockage (isHoleCell)
744 // - mark faces of cell as visible from inside/outside
745
746 //newLocationFaces = locationFaces;
747 //newHoleCell = isHoleCell;
748
749 label nChanged = 0;
750 for (const label facei : isFrontFace)
751 {
752 const label own = faceOwner[facei];
753 if (isHoleCell[own])
754 {
755 const label ownDist = surfaceCellDist[own];
756 if (ownDist >= surfaceDist)
757 {
758 //newHoleCell.unset(own);
759 isHoleCell.unset(own);
760 nChanged++;
761
762 const cell& cFaces = mesh_.cells()[own];
763 // Set corresponding bits on faces
764 const unsigned int mask = locationFaces[facei];
765 for (const label fi : cFaces)
766 {
767 //newLocationFaces[fi] |= mask;
768 locationFaces[fi] |= mask;
769 }
770 }
771 }
772 else if (mesh_.isInternalFace(facei))
773 {
774 const label nei = faceNeighbour[facei];
775 const label neiDist = surfaceCellDist[nei];
776
777 if (isHoleCell[nei] && neiDist >= surfaceDist)
778 {
779 //newHoleCell[nei] = false;
780 isHoleCell[nei] = false;
781 nChanged++;
782
783 const cell& cFaces = mesh_.cells()[nei];
784 // Set corresponding bits on faces
785 const unsigned int mask = locationFaces[facei];
786 for (const label fi : cFaces)
787 {
788 //newLocationFaces[fi] |= mask;
789 locationFaces[fi] |= mask;
790 }
791 }
792 }
793 }
794
795 reduce(nChanged, sumOp<label>());
796
797
798 if (debug)
799 {
800 Pout<< "holeToFace::findClosure :"
801 << " surfaceDist:" << surfaceDist
802 << " front:" << isFrontFace.count()
803 << " nChangedCells:" << nChanged
804 << endl;
805 }
806
807
808 if (nChanged == 0)
809 {
810 // Nothing eroded at this level. Erode cells nearer to surface.
811 --surfaceDist;
812 }
813 else
814 {
815 //locationFaces = newLocationFaces;
816 //isHoleCell = newHoleCell;
817
818 // Sync locationFaces
820 (
821 mesh_,
822 locationFaces,
823 bitOrEqOp<unsigned int>()
824 );
825
826 // Calculate new front. Never include faces that are both visible
827 // from outside and inside
828 isFrontFace = frontFaces(isSurfaceFace, locationFaces, isHoleCell);
829 }
830 }
831
832
833 // Debug: dump all the end fronts
834 //{
835 // forAll(locationCells, zonei)
836 // {
837 // const unsigned int mask = (1u << zonei);
838 //
839 // bitSet isZoneFace(mesh_.nFaces());
840 // forAll(locationFaces, facei)
841 // {
842 // if (locationFaces[facei] & mask)
843 // {
844 // isZoneFace.set(facei);
845 // }
846 // }
847 // writeFaces("isZoneFace"+Foam::name(zonei)+".obj", isZoneFace);
848 // }
849 //}
850
851
852
853 // Find faces that are connected to more than one location
854 bitSet isCommonFace(mesh_.nFaces());
855 forAll(locationFaces, facei)
856 {
857 unsigned int masks = locationFaces[facei];
858 if (masks != 0u)
859 {
860 // Count number of bits set
861 const label nSet = BitOps::bit_count(masks);
862 if (nSet >= 2)
863 {
864 isCommonFace.set(facei);
865 }
866 }
867 }
868
869 // Remove faces that are on the surface
870 for (const label facei : isCommonFace)
871 {
872 if (surfaceFaceDist[facei] == 0)
873 {
874 isCommonFace.unset(facei);
875 }
876 }
877
878 if (debug)
879 {
880 Pout<< "holeToFace::findClosure :"
881 << " closure faces:" << isCommonFace.count() << endl;
882 }
883
884 return isCommonFace;
885}
886
887
888Foam::bitSet Foam::holeToFace::erodeSet
889(
890 const bitSet& isSurfaceFace,
891 const bitSet& isSetFace
892) const
893{
894 // Detect cells with lots of faces in the set. WIP. Not parallel consistent.
895
896 const labelList& faceOwner = mesh_.faceOwner();
897 const labelList& faceNeighbour = mesh_.faceNeighbour();
898
899 bitSet isSetCell(mesh_.nCells());
900 for (const label facei : isSetFace)
901 {
902 isSetCell.set(faceOwner[facei]);
903 if (mesh_.isInternalFace(facei))
904 {
905 isSetCell.set(faceNeighbour[facei]);
906 }
907 }
908
909 // Count number of faces per cell. Decide if surface would improve by
910 // moving set
911 bitSet erodedSet(isSetFace);
912 for (const label celli : isSetCell)
913 {
914 const cell& cFaces = mesh_.cells()[celli];
915
916 label nBlockedFaces = 0;
917 label nSurfaceFaces = 0;
918 for (const label facei : cFaces)
919 {
920 if (erodedSet[facei])
921 {
922 nBlockedFaces++;
923 }
924 else if (isSurfaceFace[facei])
925 {
926 nSurfaceFaces++;
927 }
928 }
929
930 if ((nSurfaceFaces + nBlockedFaces) == cFaces.size())
931 {
932 // Single cell already disconnected by surface intersections
933 for (const label facei : cFaces)
934 {
935 erodedSet.unset(facei);
936 }
937 }
938 }
940 (
941 mesh_,
942 erodedSet,
943 andEqOp<unsigned int>()
944 );
945 //checkFaceSync("erodedSet", erodedSet);
946
947 for (const label celli : isSetCell)
948 {
949 const cell& cFaces = mesh_.cells()[celli];
950
951 label nBlockedFaces = 0;
952 for (const label facei : cFaces)
953 {
954 if (erodedSet[facei])
955 {
956 nBlockedFaces++;
957 }
958 }
959 if (nBlockedFaces >= cFaces.size()-2)
960 {
961 for (const label facei : cFaces)
962 {
963 erodedSet.flip(facei);
964 }
965 }
966 }
968 (
969 mesh_,
970 erodedSet,
971 andEqOp<unsigned int>()
972 );
973
974 if (debug)
975 {
976 Pout<< "holeToFace::erodeSet :"
977 << " starting set:" << isSetFace.count()
978 << " eroded set:" << erodedSet.count() << endl;
979 }
980
981 //checkFaceSync("erodedSet", erodedSet);
982 return erodedSet;
983}
984
985
987(
988 topoSet& set,
989 const bitSet& isSurfaceFace,
990 const bitSet& isHoleCell,
991 const bool add
992) const
993{
994 labelListList locationCells(zonePoints_.size());
995 forAll(zonePoints_, zonei)
996 {
997 const pointField& zoneLocations = zonePoints_[zonei];
998 labelList& zoneCells = locationCells[zonei];
999 zoneCells.setSize(zoneLocations.size());
1000 forAll(zoneLocations, i)
1001 {
1002 const label celli = mesh_.findCell(zoneLocations[i]);
1003 zoneCells[i] = celli;
1004
1005 // Check that cell has at least one unblocked face so front can
1006 // 'escape'.
1007 if (celli != -1)
1008 {
1009 const cell& cFaces = mesh_.cells()[celli];
1010 bool hasUnblocked = false;
1011 for (const label facei : cFaces)
1012 {
1013 if (!isSurfaceFace[facei])
1014 {
1015 hasUnblocked = true;
1016 break;
1017 }
1018 }
1019
1020 if (!hasUnblocked)
1021 {
1023 << "problem : location:" << zoneLocations[i]
1024 << " in zone:" << zonei
1025 << " is found in cell at:" << celli
1026 << mesh_.cellCentres()[celli]
1027 << " which is completely surrounded by blocked faces"
1028 << exit(FatalError);
1029 }
1030 }
1031 }
1032 }
1033
1034 bitSet isClosingFace
1035 (
1036 findClosure
1037 (
1038 isSurfaceFace, // intersected faces
1039 isHoleCell, // starting blockage
1040 locationCells // cells for zonePoints_
1041 )
1042 );
1043
1044 if (erode_)
1045 {
1046 isClosingFace = erodeSet(isSurfaceFace, isClosingFace);
1047 }
1048
1049 if (debug)
1050 {
1051 writeFaces("isClosingFace.obj", isClosingFace);
1052 //checkFaceSync("isClosingFace", isCommonFace);
1053 }
1054
1055 for (const label facei : isClosingFace)
1056 {
1057 addOrDelete(set, facei, add);
1058 }
1059}
1060
1061
1062Foam::List<Foam::pointField> Foam::holeToFace::expand(const pointField& pts)
1063{
1064 List<pointField> allPts(pts.size());
1065 forAll(pts, i)
1066 {
1067 pointField& onePt = allPts[i];
1068 onePt.setSize(1, pts[i]);
1069 }
1070 return allPts;
1071}
1072
1073
1074// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1075
1077(
1078 const polyMesh& mesh,
1079 const List<pointField>& zonePoints,
1080 const wordList& blockedFaceNames,
1081 const wordList& blockedCellNames,
1082 const bool erode
1083)
1084:
1086 zonePoints_(zonePoints),
1087 blockedFaceNames_(blockedFaceNames),
1088 blockedCellNames_(blockedCellNames),
1089 erode_(erode)
1090{}
1091
1092
1094(
1095 const polyMesh& mesh,
1096 const dictionary& dict
1097)
1098:
1100 zonePoints_(dict.get<List<pointField>>("points")),
1101 blockedFaceNames_(),
1102 blockedCellNames_(),
1103 erode_(dict.getOrDefault<bool>("erode", false))
1104{
1105 // Look for 'sets' or 'set'
1106 if (!dict.readIfPresent("faceSets", blockedFaceNames_))
1107 {
1108 blockedFaceNames_.resize(1);
1109 if (!dict.readEntry("faceSet", blockedFaceNames_.first()))
1110 {
1111 blockedFaceNames_.clear();
1112 }
1113 }
1114 if (!dict.readIfPresent("cellSets", blockedCellNames_))
1115 {
1116 blockedCellNames_.resize(1);
1117 if (!dict.readEntry("cellSet", blockedCellNames_.first()))
1118 {
1119 blockedCellNames_.clear();
1120 }
1121 }
1122}
1123
1124
1126(
1127 const polyMesh& mesh,
1128 Istream& is
1129)
1130:
1132 zonePoints_(expand(pointField(is))),
1133 blockedFaceNames_(one(), word(checkIs(is))),
1134 blockedCellNames_(),
1135 erode_(false)
1136{}
1137
1138
1139// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1140
1142(
1143 const topoSetSource::setAction action,
1144 topoSet& set
1145) const
1146{
1147 // Set of additional blocked (internal or coupled) faces
1148 bitSet isBlockedFace(mesh_.nFaces());
1149 for (const word& setName : blockedFaceNames_)
1150 {
1151 const faceSet loadedSet(mesh_, setName);
1152 isBlockedFace.set(loadedSet.toc());
1153 }
1154
1155 // Optional initial blocked cells
1156 bitSet isCandidateCell(mesh_.nCells());
1157 if (blockedFaceNames_.size())
1158 {
1159 for (const word& setName : blockedCellNames_)
1160 {
1161 const cellSet loadedSet(mesh_, setName);
1162 isCandidateCell.set(loadedSet.toc());
1163 }
1164 }
1165 else
1166 {
1167 isCandidateCell = true;
1168 }
1169
1170 if (action == topoSetSource::ADD || action == topoSetSource::NEW)
1171 {
1172 if (verbose_)
1173 {
1174 Info<< " Adding all faces to disconnect regions "
1175 << flatOutput(zonePoints_) << " ..." << endl;
1176 }
1177
1178 combine(set, isBlockedFace, isCandidateCell, true);
1179 }
1180 else if (action == topoSetSource::SUBTRACT)
1181 {
1182 if (verbose_)
1183 {
1184 Info<< " Removing all faces to disconnect regions "
1185 << flatOutput(zonePoints_) << " ..." << endl;
1186 }
1187
1188 combine(set, isBlockedFace, isCandidateCell, false);
1189 }
1190}
1191
1192
1194(
1195 const polyMesh& mesh,
1196 const List<pointField>& zonePoints,
1197 const labelList& blockedFaces,
1198 const globalIndex& globalBlockedFaces,
1199 const bool erode,
1200
1201 labelList& closureFaces, // local faces to close gap
1202 labelList& closureToBlocked
1203)
1204{
1205 if (blockedFaces.size() != globalBlockedFaces.localSize())
1206 {
1207 FatalErrorInFunction << "problem : blockedFaces:" << blockedFaces.size()
1208 << " globalBlockedFaces:" << globalBlockedFaces.localSize()
1209 << exit(FatalError);
1210 }
1211
1212
1213 // Calculate faces needed to close hole (closureFaces)
1214 {
1215 const holeToFace faceSetSource
1216 (
1217 mesh,
1218 zonePoints,
1219 wordList(0),
1220 wordList(0),
1221 erode //false
1222 );
1223 faceSet closureFaceSet(mesh, "calcClosure", 256);
1224
1225 const bitSet isBlockedFace(mesh.nFaces(), blockedFaces);
1226 const bitSet isActiveCell(mesh.nCells(), true);
1227
1228 faceSetSource.combine
1229 (
1230 closureFaceSet,
1231 isBlockedFace,
1232 isActiveCell,
1233 true
1234 );
1235
1236 closureFaces = closureFaceSet.sortedToc();
1237 }
1238
1239
1240 if (returnReduce(closureFaces.size(), sumOp<label>()) == 0)
1241 {
1242 closureToBlocked.clear();
1243 return autoPtr<mapDistribute>(nullptr);
1244 }
1245
1246
1247 //- Seed edges of closureFaces patch with (global) index of blockedFace
1248
1249 const indirectPrimitivePatch pp
1250 (
1251 IndirectList<face>(mesh.faces(), closureFaces),
1252 mesh.points()
1253 );
1254 const edgeList& edges = pp.edges();
1255 const labelList& mp = pp.meshPoints();
1256 const label nBndEdges = pp.nEdges() - pp.nInternalEdges();
1257
1258 // For all faces in blockedFaces mark the edge with a face. No special
1259 // handling for multiple faces sharing the edge - first one wins
1260 EdgeMap<label> edgeMap(pp.nEdges());
1261 forAll(blockedFaces, i)
1262 {
1263 const label globalBlockedi = globalBlockedFaces.toGlobal(i);
1264 const label facei = blockedFaces[i];
1265 const face& f = mesh.faces()[facei];
1266 forAll(f, fp)
1267 {
1268 label nextFp = f.fcIndex(fp);
1269 edgeMap.insert(edge(f[fp], f[nextFp]), globalBlockedi);
1270 }
1271 }
1273
1274
1275
1276 // Seed
1277 DynamicList<label> initialEdges(2*nBndEdges);
1279 initialEdgesInfo(2*nBndEdges);
1280 forAll(edges, edgei)
1281 {
1282 const edge& e = edges[edgei];
1283 const edge meshE = edge(mp[e[0]], mp[e[1]]);
1284
1285 auto iter = edgeMap.cfind(meshE);
1286 if (iter.found())
1287 {
1288 // Found edge on patch connected to blocked face. Seed with the
1289 // (global) index of that blocked face
1290
1291 initialEdges.append(edgei);
1292 initialEdgesInfo.append
1293 (
1295 (
1296 0, // distance
1297 iter() // globalBlockedi
1298 )
1299 );
1300 }
1301 }
1302
1303 // Data on all edges and faces
1305 (
1306 pp.nEdges()
1307 );
1309 (
1310 pp.size()
1311 );
1312
1313 // Walk
1315 <
1318 > calc
1319 (
1320 mesh,
1321 pp,
1322 initialEdges,
1323 initialEdgesInfo,
1324 allEdgeInfo,
1325 allFaceInfo,
1327 );
1328
1329
1330 // Per closure face the seed face
1331 closureToBlocked.setSize(pp.size());
1332 closureToBlocked = -1;
1333 forAll(allFaceInfo, facei)
1334 {
1335 if (allFaceInfo[facei].valid(calc.data()))
1336 {
1337 closureToBlocked[facei] = allFaceInfo[facei].data();
1338 }
1339 }
1340 // Above wave only guarantees unique data on coupled edges, not on
1341 // coupled faces (?) so explicitly sync faces
1342 {
1343 labelList syncFld(mesh.nFaces(), -1);
1344 UIndirectList<label>(syncFld, pp.addressing()) = closureToBlocked;
1346 closureToBlocked = UIndirectList<label>(syncFld, pp.addressing());
1347 }
1348
1349 List<Map<label>> compactMap;
1351 (
1352 globalBlockedFaces,
1353 closureToBlocked,
1354 compactMap
1355 );
1356}
1357
1358
1359// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
Map from edge (expressed as its endpoints) to value. For easier forward declaration it is currently i...
Definition: EdgeMap.H:54
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:122
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition: HashTableI.H:141
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
A List with indirect addressing.
Definition: IndirectList.H:119
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Wave propagation of information along patch. Every iteration information goes through one layer of fa...
A list of faces which address into the list of points.
label nEdges() const
Number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
fileName timePath() const
Return current time path.
Definition: Time.H:375
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
T & first()
Return the first element of the list.
Definition: UListI.H:202
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:237
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:590
A collection of cell labels.
Definition: cellSet.H:54
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
For use with PatchEdgeFaceWave. Determines topological distance to starting edges....
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
A list of face labels.
Definition: faceSet.H:54
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
label localSize() const
My local size.
Definition: globalIndexI.H:207
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:260
A topoSetFaceSource to select a set of faces that closes a hole i.e. disconnects zones (specified by ...
Definition: holeToFace.H:95
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
Apply specified action to the topoSet.
Definition: holeToFace.C:1142
static autoPtr< mapDistribute > calcClosure(const polyMesh &mesh, const List< pointField > &zonePoints, const labelList &blockedFaces, const globalIndex &globalBlockedFaces, const bool erode, labelList &closureFaces, labelList &closureToBlocked)
Optional direct use to generate the set of faces and the method to.
Definition: holeToFace.C:1194
void combine(topoSet &set, const bitSet &isBlockedFace, const bitSet &isActiveCell, const bool add) const
Optional direct use to generate a faceSet.
Definition: holeToFace.C:987
const Time & time() const noexcept
Return time registry.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
static void syncEdgeMap(const polyMesh &mesh, EdgeMap< T > &edgeValues, const CombineOp &cop, const TransformOp &top)
Synchronize values on selected edges.
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:396
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:445
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
The topoSetFaceSource is a intermediate class for handling topoSet sources for selecting faces.
Class with constructor to add usage string to table.
Base class of a source for a topoSet.
Definition: topoSetSource.H:68
setAction
Enumeration defining various actions.
@ SUBTRACT
Subtract elements from current set.
@ ADD
Add elements to current set.
@ NEW
Create a new set and ADD elements to it.
const polyMesh & mesh_
Reference to the mesh.
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:67
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
bool
Definition: EEqn.H:20
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nBndEdges(UPstream::listGatherValues< label >(aMesh.nBoundaryEdges() - nLocalProcEdges))
unsigned int bit_count(UIntType x)
Count arbitrary number of bits (of an integral type)
Definition: BitOps.H:171
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
List< label > labelList
A List of labels.
Definition: List.H:66
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:515
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
3D tensor transformation operations.