fvMeshSubset.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) 2015-2022 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
29#include "fvMeshSubset.H"
30#include "BitOps.H"
31#include "Pstream.H"
32#include "cyclicPolyPatch.H"
33#include "emptyPolyPatch.H"
34#include "processorPolyPatch.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
39
40
41// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
42
43namespace
44{
45
46// Mark faces/points (with 0) in labelList
47inline void markUsed
48(
49 const Foam::labelUList& locations,
51)
52{
53 for (auto idx : locations)
54 {
55 map[idx] = 0;
56 }
57}
58
59} // End anonymous namespace
60
61
62namespace Foam
63{
64
65// Perform a subset of a subset
67(
68 const label nElems,
69 const labelUList& selectedElements, // First subset
70 const labelUList& subsetMap // Subset within first subset
71)
72{
73 if (selectedElements.empty() || subsetMap.empty())
74 {
75 // Trivial case
76 return labelList();
77 }
78
79 // Mark selected elements.
80 const bitSet selected(nElems, selectedElements);
81
82 // Count subset of selected elements
83 label n = 0;
84 forAll(subsetMap, i)
85 {
86 if (selected[subsetMap[i]])
87 {
88 ++n;
89 }
90 }
91
92 // Collect selected elements
93 labelList subsettedElements(n);
94 n = 0;
95
96 forAll(subsetMap, i)
97 {
98 if (selected[subsetMap[i]])
99 {
100 subsettedElements[n] = i;
101 ++n;
102 }
103 }
104
105 return subsettedElements;
106}
107
108} // End namespace Foam
109
110
111// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
112
114{
115 if (!subMeshPtr_)
116 {
118 << "Mesh is not subsetted!" << nl
119 << abort(FatalError);
120
121 return false;
122 }
123
124 return true;
125}
126
127
128void Foam::fvMeshSubset::calcFaceFlipMap() const
129{
130 const labelList& subToBaseFace = faceMap();
131 const labelList& subToBaseCell = cellMap();
132
133 faceFlipMapPtr_.reset(new labelList(subToBaseFace.size()));
134 auto& faceFlipMap = *faceFlipMapPtr_;
135
136 // Only exposed internal faces might be flipped (since we don't do
137 // any cell renumbering, just compacting)
138 const label subInt = subMesh().nInternalFaces();
139
140 const labelList& subOwn = subMesh().faceOwner();
141 const labelList& own = baseMesh_.faceOwner();
142
143 for (label subFaceI = 0; subFaceI < subInt; ++subFaceI)
144 {
145 faceFlipMap[subFaceI] = subToBaseFace[subFaceI]+1;
147 for (label subFaceI = subInt; subFaceI < subOwn.size(); ++subFaceI)
148 {
149 const label faceI = subToBaseFace[subFaceI];
150 if (subToBaseCell[subOwn[subFaceI]] == own[faceI])
151 {
152 faceFlipMap[subFaceI] = faceI+1;
154 else
155 {
156 faceFlipMap[subFaceI] = -faceI-1;
157 }
158 }
159}
160
161
162void Foam::fvMeshSubset::doCoupledPatches
164 const bool syncPar,
165 labelUList& nCellsUsingFace
166) const
167{
168 // Synchronize facesToSubset on both sides of coupled patches.
169 // Marks faces that become 'uncoupled' with 3.
170
171 const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
172
173 label nUncoupled = 0;
174
175 if (syncPar && Pstream::parRun())
176 {
178
179 // Send face usage across processor patches
180 for (const polyPatch& pp : oldPatches)
181 {
182 const auto* procPatch = isA<processorPolyPatch>(pp);
184 if (procPatch)
185 {
186 const label nbrProci = procPatch->neighbProcNo();
187
188 UOPstream toNeighbour(nbrProci, pBufs);
189
190 if (!nCellsUsingFace.empty())
191 {
192 toNeighbour <<
193 SubList<label>(nCellsUsingFace, pp.size(), pp.start());
194 }
195 else
196 {
197 toNeighbour << labelList();
198 }
199 }
200 }
201
202 pBufs.finishedSends();
203
204 // Receive face usage count and check for faces that become uncoupled.
205 for (const polyPatch& pp : oldPatches)
206 {
207 const auto* procPatch = isA<processorPolyPatch>(pp);
208
209 if (procPatch)
210 {
211 const label nbrProci = procPatch->neighbProcNo();
212
213 UIPstream fromNeighbour(nbrProci, pBufs);
214
215 const labelList nbrList(fromNeighbour);
216
217 // Combine with this side.
218
219 if (!nCellsUsingFace.empty())
220 {
221 const labelList& nbrCellsUsingFace(nbrList);
222
223 // Combine with this side.
224
225 forAll(pp, i)
226 {
227 if
228 (
229 nCellsUsingFace[pp.start()+i] == 1
230 && nbrCellsUsingFace[i] == 0
231 )
232 {
233 // Face's neighbour is no longer there.
234 // Mark face off as coupled
235 nCellsUsingFace[pp.start()+i] = 3;
236 ++nUncoupled;
237 }
238 }
239 }
240 }
241 }
242 }
243
244 // Do same for cyclics.
245 for (const polyPatch& pp : oldPatches)
246 {
247 const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(pp);
248
249 if (cpp && !nCellsUsingFace.empty())
250 {
251 const auto& cycPatch = *cpp;
252
253 forAll(cycPatch, i)
254 {
255 label thisFacei = cycPatch.start() + i;
256 label otherFacei = cycPatch.transformGlobalFace(thisFacei);
257
258 if
259 (
260 nCellsUsingFace[thisFacei] == 1
261 && nCellsUsingFace[otherFacei] == 0
262 )
263 {
264 nCellsUsingFace[thisFacei] = 3;
265 ++nUncoupled;
266 }
267 }
268 }
269 }
270
271 if (syncPar)
272 {
273 reduce(nUncoupled, sumOp<label>());
274 }
275
276 if (nUncoupled)
277 {
278 Info<< "Uncoupled " << nUncoupled << " faces on coupled patches. "
279 << "(processorPolyPatch, cyclicPolyPatch)" << endl;
280 }
281}
282
283
284void Foam::fvMeshSubset::subsetZones()
285{
286 // Keep all zones, even if zero size.
287
288 #ifdef FULLDEBUG
289 checkHasSubMesh();
290 #endif
291
292 auto& newSubMesh = subMeshPtr_();
293
294 // PointZones
295
296 const pointZoneMesh& pointZones = baseMesh().pointZones();
297
298 List<pointZone*> pZones(pointZones.size());
299
300 forAll(pointZones, zonei)
301 {
302 const pointZone& pz = pointZones[zonei];
303
304 pZones[zonei] = new pointZone
305 (
306 pz.name(),
307 subsetSubset(baseMesh().nPoints(), pz, pointMap()),
308 zonei,
309 newSubMesh.pointZones()
310 );
311 }
312
313
314 // FaceZones
315 // Do we need to remove zones where the side we're interested in
316 // no longer exists? Guess not.
317
318 const faceZoneMesh& faceZones = baseMesh().faceZones();
319
320 List<faceZone*> fZones(faceZones.size());
321
322 forAll(faceZones, zonei)
323 {
324 const faceZone& fz = faceZones[zonei];
325
326 // Expand faceZone to full mesh
327 // +1 : part of faceZone, flipped
328 // -1 : ,, , unflipped
329 // 0 : not part of faceZone
330 labelList zone(baseMesh().nFaces(), Zero);
331 forAll(fz, j)
332 {
333 zone[fz[j]] = (fz.flipMap()[j] ? 1 : -1);
334 }
335
336 // Select faces
337 label nSub = 0;
338 forAll(faceMap(), j)
339 {
340 if (zone[faceMap()[j]] != 0)
341 {
342 ++nSub;
343 }
344 }
345 labelList subAddressing(nSub);
346 boolList subFlipStatus(nSub);
347 nSub = 0;
348 forAll(faceMap(), subFacei)
349 {
350 const label meshFacei = faceMap()[subFacei];
351 if (zone[meshFacei] != 0)
352 {
353 subAddressing[nSub] = subFacei;
354 const label subOwner = subMesh().faceOwner()[subFacei];
355 const label baseOwner = baseMesh().faceOwner()[meshFacei];
356 // If subowner is the same cell as the base keep the flip status
357 const bool sameOwner = (cellMap()[subOwner] == baseOwner);
358 const bool flip = (zone[meshFacei] == 1);
359 subFlipStatus[nSub] = (sameOwner == flip);
360
361 ++nSub;
362 }
363 }
364
365 fZones[zonei] = new faceZone
366 (
367 fz.name(),
368 subAddressing,
369 subFlipStatus,
370 zonei,
371 newSubMesh.faceZones()
372 );
373 }
374
375
376 // Cell Zones
377
378 const cellZoneMesh& cellZones = baseMesh().cellZones();
379
380 List<cellZone*> cZones(cellZones.size());
381
382 forAll(cellZones, zonei)
383 {
384 const cellZone& cz = cellZones[zonei];
385
386 cZones[zonei] = new cellZone
387 (
388 cz.name(),
389 subsetSubset(baseMesh().nCells(), cz, cellMap()),
390 zonei,
391 newSubMesh.cellZones()
392 );
393 }
394
395 // Add the zones
396 newSubMesh.addZones(pZones, fZones, cZones);
397}
398
399
400// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
401
403:
404 baseMesh_(baseMesh),
405 subMeshPtr_(nullptr),
406 faceFlipMapPtr_(nullptr),
407 pointMap_(),
408 faceMap_(),
409 cellMap_(),
410 patchMap_()
411{}
412
413
415:
416 fvMeshSubset(baseMesh)
417{
418 reset(Foam::zero{});
419}
420
421
423(
424 const fvMesh& baseMesh,
425 const bitSet& selectedCells,
426 const label patchID,
427 const bool syncPar
428)
429:
430 fvMeshSubset(baseMesh)
431{
432 reset(selectedCells, patchID, syncPar);
433}
434
435
437(
438 const fvMesh& baseMesh,
439 const labelUList& selectedCells,
440 const label patchID,
441 const bool syncPar
442)
443:
444 fvMeshSubset(baseMesh)
445{
446 reset(selectedCells, patchID, syncPar);
447}
448
449
451(
452 const fvMesh& baseMesh,
453 const labelHashSet& selectedCells,
454 const label patchID,
455 const bool syncPar
456)
457:
458 fvMeshSubset(baseMesh)
459{
460 reset(selectedCells, patchID, syncPar);
461}
462
463
465(
466 const fvMesh& baseMesh,
467 const label regioni,
468 const labelUList& regions,
469 const label patchID,
470 const bool syncPar
471)
472:
473 fvMeshSubset(baseMesh)
474{
475 reset(regioni, regions, patchID, syncPar);
476}
477
478
479// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
480
482{
483 subMeshPtr_.reset(nullptr);
484 faceFlipMapPtr_.reset(nullptr);
485
486 pointMap_.clear();
487 faceMap_.clear();
488 cellMap_.clear();
489 patchMap_.clear();
490}
491
492
494{
495 clear();
496}
497
498
500(
501 autoPtr<fvMesh>&& subMeshPtr,
502 labelList&& pointMap,
504 labelList&& cellMap,
505 labelList&& patchMap
506)
507{
508 subMeshPtr_.reset(std::move(subMeshPtr));
509 faceFlipMapPtr_.reset(nullptr);
510
511 pointMap_ = std::move(pointMap);
512 faceMap_ = std::move(faceMap);
513 cellMap_ = std::move(cellMap);
514 patchMap_ = std::move(patchMap);
515
516 // Sanity
517 if (!subMeshPtr_)
518 {
519 clear();
520 }
521}
522
523
525{
526 clear();
527
528 // Create zero-sized subMesh
529 subMeshPtr_.reset
530 (
531 new fvMesh
532 (
534 (
535 baseMesh_.name(),
536 baseMesh_.time().timeName(),
537 baseMesh_.time(),
538 IOobject::NO_READ, // Do not read any dictionaries
540 ),
541 baseMesh_, // Get dictionaries from base mesh
542 Foam::zero{} // zero-sized
543 // Uses syncPar (bounds) - should generally be OK
544 )
545 );
546 auto& newSubMesh = subMeshPtr_();
547
548
549 // Clone non-processor patches
550 {
551 const polyBoundaryMesh& oldBoundary = baseMesh_.boundaryMesh();
552 const polyBoundaryMesh& newBoundary = newSubMesh.boundaryMesh();
553
554 polyPatchList newPatches(oldBoundary.nNonProcessor());
555
556 patchMap_ = identity(newPatches.size());
557
558 forAll(newPatches, patchi)
559 {
560 newPatches.set
561 (
562 patchi,
563 oldBoundary[patchi].clone
564 (
565 newBoundary,
566 patchi,
567 0, // patch size
568 0 // patch start
569 )
570 );
571 }
572
573 newSubMesh.addFvPatches(newPatches);
574 }
575
576
577 // Add the zones
578 subsetZones();
579}
580
581
583(
584 const bitSet& selectedCells,
585 const label patchID,
586 const bool syncPar
587)
588{
589 // Clear all old maps and pointers
590 clear();
591
592 const cellList& oldCells = baseMesh().cells();
593 const faceList& oldFaces = baseMesh().faces();
594 const pointField& oldPoints = baseMesh().points();
595 const labelList& oldOwner = baseMesh().faceOwner();
596 const labelList& oldNeighbour = baseMesh().faceNeighbour();
597 const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
598 const label oldNInternalFaces = baseMesh().nInternalFaces();
599
600 // Initial check on patches before doing anything time consuming.
601
602 label wantedPatchID = patchID;
603
604 if (wantedPatchID == -1)
605 {
606 // No explicit patch specified. Put in oldInternalFaces patch.
607 // Check if patch with this name already exists.
608 wantedPatchID = oldPatches.findPatchID(exposedPatchName);
609 }
610 else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
611 {
613 << "Non-existing patch index " << wantedPatchID << endl
614 << "Should be between 0 and " << oldPatches.size()-1
615 << abort(FatalError);
616 }
617
618 // The selected cells - sorted in ascending order
619 cellMap_ = selectedCells.sortedToc();
620
621 // The selectedCells should normally be in the [0,nCells) range,
622 // but it is better not to trust that.
623 {
624 label len = cellMap_.size();
625 for
626 (
627 label i = len-1;
628 i >= 0 && (cellMap_[i] >= oldCells.size());
629 --i
630 )
631 {
632 len = i;
633 }
634 cellMap_.resize(len);
635 }
636
637
638 // Mark all used faces. Count number of cells using them
639 // 0: face not used anymore
640 // 1: face used by one cell, face becomes/stays boundary face
641 // 2: face still used and remains internal face
642 // 3: face coupled and used by one cell only (so should become normal,
643 // non-coupled patch face)
644 //
645 // Note that this is not really necessary - but means we can size things
646 // correctly. Also makes handling coupled faces much easier.
647
648 labelList nCellsUsingFace(oldFaces.size(), Zero);
649
650 label nFacesInSet = 0;
651
652 forAll(oldFaces, oldFacei)
653 {
654 bool faceUsed = false;
655
656 if (selectedCells.test(oldOwner[oldFacei]))
657 {
658 ++nCellsUsingFace[oldFacei];
659 faceUsed = true;
660 }
661
662 if
663 (
664 baseMesh().isInternalFace(oldFacei)
665 && selectedCells.test(oldNeighbour[oldFacei])
666 )
667 {
668 ++nCellsUsingFace[oldFacei];
669 faceUsed = true;
670 }
671
672 if (faceUsed)
673 {
674 ++nFacesInSet;
675 }
676 }
677 faceMap_.resize(nFacesInSet);
678
679 // Handle coupled faces. Modifies patch faces to be uncoupled to 3.
680 doCoupledPatches(syncPar, nCellsUsingFace);
681
682
683 // See which patch to use for exposed internal faces.
684 label oldInternalPatchID = 0;
685
686 // Insert faces before which patch
687 label nextPatchID = oldPatches.size();
688
689 // old to new patches
690 labelList globalPatchMap(oldPatches.size());
691
692 // New patch size
693 label nbSize = oldPatches.size();
694
695 if (wantedPatchID == -1)
696 {
697 // Create 'oldInternalFaces' patch at the end (or before
698 // processorPatches)
699 // and put all exposed internal faces in there.
700
701 forAll(oldPatches, patchi)
702 {
703 if (isA<processorPolyPatch>(oldPatches[patchi]))
704 {
705 nextPatchID = patchi;
706 break;
707 }
708 ++oldInternalPatchID;
709 }
710
711 ++nbSize;
712
713 // adapt old to new patches for inserted patch
714 for (label oldPatchi = 0; oldPatchi < nextPatchID; oldPatchi++)
715 {
716 globalPatchMap[oldPatchi] = oldPatchi;
717 }
718 for
719 (
720 label oldPatchi = nextPatchID;
721 oldPatchi < oldPatches.size();
722 oldPatchi++
723 )
724 {
725 globalPatchMap[oldPatchi] = oldPatchi+1;
726 }
727 }
728 else
729 {
730 oldInternalPatchID = wantedPatchID;
731 nextPatchID = wantedPatchID+1;
732
733 // old to new patches
734 globalPatchMap = identity(oldPatches.size());
735 }
736
737 labelList boundaryPatchSizes(nbSize, Zero);
738
739
740 // Make a global-to-local point map
741 labelList globalPointMap(oldPoints.size(), -1);
742 labelList globalFaceMap(oldFaces.size(), -1);
743
744 label facei = 0;
745
746 // 1. Pick up all preserved internal faces.
747 for (label oldFacei = 0; oldFacei < oldNInternalFaces; ++oldFacei)
748 {
749 if (nCellsUsingFace[oldFacei] == 2)
750 {
751 globalFaceMap[oldFacei] = facei;
752 faceMap_[facei++] = oldFacei;
753
754 // Mark all points from the face
755 markUsed(oldFaces[oldFacei], globalPointMap);
756 }
757 }
758
759 // These are all the internal faces in the mesh.
760 const label nInternalFaces = facei;
761
762 // 2. Boundary faces up to where we want to insert old internal faces
763 for
764 (
765 label oldPatchi = 0;
766 oldPatchi < oldPatches.size()
767 && oldPatchi < nextPatchID;
768 oldPatchi++
769 )
770 {
771 const polyPatch& oldPatch = oldPatches[oldPatchi];
772
773 label oldFacei = oldPatch.start();
774
775 forAll(oldPatch, i)
776 {
777 if (nCellsUsingFace[oldFacei] == 1)
778 {
779 // Boundary face is kept.
780
781 // Mark face and increment number of points in set
782 globalFaceMap[oldFacei] = facei;
783 faceMap_[facei++] = oldFacei;
784
785 // Mark all points from the face
786 markUsed(oldFaces[oldFacei], globalPointMap);
787
788 // Increment number of patch faces
789 ++boundaryPatchSizes[globalPatchMap[oldPatchi]];
790 }
791 ++oldFacei;
792 }
793 }
794
795 // 3a. old internal faces that have become exposed.
796 for (label oldFacei = 0; oldFacei < oldNInternalFaces; ++oldFacei)
797 {
798 if (nCellsUsingFace[oldFacei] == 1)
799 {
800 globalFaceMap[oldFacei] = facei;
801 faceMap_[facei++] = oldFacei;
802
803 // Mark all points from the face
804 markUsed(oldFaces[oldFacei], globalPointMap);
805
806 // Increment number of patch faces
807 ++boundaryPatchSizes[oldInternalPatchID];
808 }
809 }
810
811 // 3b. coupled patch faces that have become uncoupled.
812 for
813 (
814 label oldFacei = oldNInternalFaces;
815 oldFacei < oldFaces.size();
816 oldFacei++
817 )
818 {
819 if (nCellsUsingFace[oldFacei] == 3)
820 {
821 globalFaceMap[oldFacei] = facei;
822 faceMap_[facei++] = oldFacei;
823
824 // Mark all points from the face
825 markUsed(oldFaces[oldFacei], globalPointMap);
826
827 // Increment number of patch faces
828 ++boundaryPatchSizes[oldInternalPatchID];
829 }
830 }
831
832 // 4. Remaining boundary faces
833 for
834 (
835 label oldPatchi = nextPatchID;
836 oldPatchi < oldPatches.size();
837 oldPatchi++
838 )
839 {
840 const polyPatch& oldPatch = oldPatches[oldPatchi];
841
842 label oldFacei = oldPatch.start();
843
844 forAll(oldPatch, i)
845 {
846 if (nCellsUsingFace[oldFacei] == 1)
847 {
848 // Boundary face is kept.
849
850 // Mark face and increment number of points in set
851 globalFaceMap[oldFacei] = facei;
852 faceMap_[facei++] = oldFacei;
853
854 // Mark all points from the face
855 markUsed(oldFaces[oldFacei], globalPointMap);
856
857 // Increment number of patch faces
858 ++boundaryPatchSizes[globalPatchMap[oldPatchi]];
859 }
860 ++oldFacei;
861 }
862 }
863
864 if (facei != nFacesInSet)
865 {
867 << "Problem" << abort(FatalError);
868 }
869
870
871 // Grab the points map
872 label nPointsInSet = 0;
873
874 forAll(globalPointMap, pointi)
875 {
876 if (globalPointMap[pointi] != -1)
877 {
878 ++nPointsInSet;
879 }
880 }
881 pointMap_.setSize(nPointsInSet);
882
883 nPointsInSet = 0;
884
885 forAll(globalPointMap, pointi)
886 {
887 if (globalPointMap[pointi] != -1)
888 {
889 pointMap_[nPointsInSet] = pointi;
890 globalPointMap[pointi] = nPointsInSet;
891 ++nPointsInSet;
892 }
893 }
894
895
896 //Pout<< "Number of points,faces,cells in new mesh : "
897 // << pointMap_.size() << ' '
898 // << faceMap_.size() << ' '
899 // << cellMap_.size() << nl;
900
901
902 // Make a new mesh
903
904 //
905 // Create new points
906 //
907 pointField newPoints(oldPoints, pointMap_);
908
909
910 //
911 // Create new faces
912 //
913
914 faceList newFaces(faceMap_.size());
915 {
916 auto iter = newFaces.begin();
917 const auto& renumbering = globalPointMap;
918
919 // Internal faces
920 for (label facei = 0; facei < nInternalFaces; ++facei)
921 {
922 face& newItem = *iter;
923 ++iter;
924
925 const face& oldItem = oldFaces[faceMap_[facei]];
926
927 // Copy relabelled
928 newItem.resize(oldItem.size());
929
930 forAll(oldItem, i)
931 {
932 newItem[i] = renumbering[oldItem[i]];
933 }
934 }
935
936 // Boundary faces - may need to be flipped
937 for (label facei = nInternalFaces; facei < faceMap_.size(); ++facei)
938 {
939 const label oldFacei = faceMap_[facei];
940
941 face& newItem = *iter;
942 ++iter;
943
944 const face oldItem =
945 (
946 (
947 baseMesh().isInternalFace(oldFacei)
948 && selectedCells.test(oldNeighbour[oldFacei])
949 && !selectedCells.test(oldOwner[oldFacei])
950 )
951 // Face flipped to point outwards
952 ? oldFaces[oldFacei].reverseFace()
953 : oldFaces[oldFacei]
954 );
955
956 // Copy relabelled
957 newItem.resize(oldItem.size());
958
959 forAll(oldItem, i)
960 {
961 newItem[i] = renumbering[oldItem[i]];
962 }
963 }
964 }
965
966
967 //
968 // Create new cells
969 //
970
971 cellList newCells(cellMap_.size());
972 {
973 auto iter = newCells.begin();
974 const auto& renumbering = globalFaceMap;
975
976 for (const label oldCelli : cellMap_)
977 {
978 cell& newItem = *iter;
979 ++iter;
980
981 const labelList& oldItem = oldCells[oldCelli];
982
983 // Copy relabelled
984 newItem.resize(oldItem.size());
985
986 forAll(oldItem, i)
987 {
988 newItem[i] = renumbering[oldItem[i]];
989 }
990 }
991 }
992
993
994 // Make a new mesh
995 //
996 // Note that mesh gets registered with same name as original mesh.
997 // This is not proper but cannot be avoided since otherwise
998 // surfaceInterpolation cannot find its fvSchemes.
999 // It will try to read, for example. "system/region0SubSet/fvSchemes"
1000 //
1001 subMeshPtr_ = autoPtr<fvMesh>::New
1002 (
1003 IOobject
1004 (
1005 baseMesh_.name(),
1006 baseMesh_.time().timeName(),
1007 baseMesh_.time(),
1008 IOobject::NO_READ, // Do not read any dictionaries
1010 ),
1011 baseMesh_, // Get dictionaries from base mesh
1012 std::move(newPoints),
1013 std::move(newFaces),
1014 std::move(newCells),
1015 syncPar // parallel synchronisation
1016 );
1017
1018
1019 // Add old patches
1020 polyPatchList newBoundary(nbSize);
1021 patchMap_.resize(nbSize);
1022 label nNewPatches = 0;
1023 label patchStart = nInternalFaces;
1024
1025
1026 // For parallel: only remove patch if none of the processors has it.
1027 // This only gets done for patches before the one being inserted
1028 // (so patches < nextPatchID)
1029
1030 // Get sum of patch sizes. Zero if patch can be deleted.
1031 labelList globalPatchSizes(boundaryPatchSizes);
1032 globalPatchSizes.setSize(nextPatchID);
1033
1034 if (syncPar && Pstream::parRun())
1035 {
1036 // Get patch names (up to nextPatchID)
1038 patchNames[Pstream::myProcNo()] = oldPatches.names();
1039 patchNames[Pstream::myProcNo()].setSize(nextPatchID);
1041
1042 // Get patch sizes (up to nextPatchID).
1043 // Note that up to nextPatchID the globalPatchMap is an identity so
1044 // no need to index through that.
1046
1047 // Now all processors have all the patchnames.
1048 // Decide: if all processors have the same patch names and size is zero
1049 // everywhere remove the patch.
1050 bool samePatches = true;
1051
1052 for (label proci = 1; proci < patchNames.size(); ++proci)
1053 {
1054 if (patchNames[proci] != patchNames[0])
1055 {
1056 samePatches = false;
1057 break;
1058 }
1059 }
1060
1061 if (!samePatches)
1062 {
1063 // Patchnames not sync on all processors so disable removal of
1064 // zero sized patches.
1065 globalPatchSizes = labelMax;
1066 }
1067 }
1068
1069
1070 // Old patches
1071
1072 for
1073 (
1074 label oldPatchi = 0;
1075 oldPatchi < oldPatches.size()
1076 && oldPatchi < nextPatchID;
1077 oldPatchi++
1078 )
1079 {
1080 const label newSize = boundaryPatchSizes[globalPatchMap[oldPatchi]];
1081
1082 if (oldInternalPatchID != oldPatchi)
1083 {
1084 // Pure subset of patch faces (no internal faces added to this
1085 // patch). Can use mapping.
1086 labelList map(newSize);
1087 for (label patchFacei = 0; patchFacei < newSize; patchFacei++)
1088 {
1089 const label facei = patchStart+patchFacei;
1090 const label oldFacei = faceMap_[facei];
1091 map[patchFacei] = oldPatches[oldPatchi].whichFace(oldFacei);
1092 }
1093
1094 newBoundary.set
1095 (
1096 nNewPatches,
1097 oldPatches[oldPatchi].clone
1098 (
1099 subMeshPtr_().boundaryMesh(),
1100 nNewPatches,
1101 map,
1102 patchStart
1103 )
1104 );
1105 }
1106 else
1107 {
1108 // Clone (even if 0 size)
1109 newBoundary.set
1110 (
1111 nNewPatches,
1112 oldPatches[oldPatchi].clone
1113 (
1114 subMeshPtr_().boundaryMesh(),
1115 nNewPatches,
1116 newSize,
1117 patchStart
1118 )
1119 );
1120 }
1121
1122 patchStart += newSize;
1123 patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1124 ++nNewPatches;
1125 }
1126
1127 // Inserted patch
1128
1129 if (wantedPatchID == -1)
1130 {
1131 label oldInternalSize = boundaryPatchSizes[oldInternalPatchID];
1132
1133 if (syncPar)
1134 {
1135 reduce(oldInternalSize, sumOp<label>());
1136 }
1137
1138 // Newly created patch so is at end. Check if any faces in it.
1139 if (oldInternalSize > 0)
1140 {
1141 newBoundary.set
1142 (
1143 nNewPatches,
1144 new emptyPolyPatch
1145 (
1146 exposedPatchName,
1147 boundaryPatchSizes[oldInternalPatchID],
1148 patchStart,
1149 nNewPatches,
1150 subMeshPtr_().boundaryMesh(),
1152 )
1153 );
1154
1155 //Pout<< " " << exposedPatchName << " : "
1156 // << boundaryPatchSizes[oldInternalPatchID] << endl;
1157
1158 // The index for the first patch is -1 as it originates from
1159 // the internal faces
1160 patchStart += boundaryPatchSizes[oldInternalPatchID];
1161 patchMap_[nNewPatches] = -1;
1162 ++nNewPatches;
1163 }
1164 }
1165
1166 // Old patches
1167
1168 for
1169 (
1170 label oldPatchi = nextPatchID;
1171 oldPatchi < oldPatches.size();
1172 oldPatchi++
1173 )
1174 {
1175 const label newSize = boundaryPatchSizes[globalPatchMap[oldPatchi]];
1176
1177 if (oldInternalPatchID != oldPatchi)
1178 {
1179 // Pure subset of patch faces (no internal faces added to this
1180 // patch). Can use mapping.
1181 labelList map(newSize);
1182 for (label patchFacei = 0; patchFacei < newSize; patchFacei++)
1183 {
1184 const label facei = patchStart+patchFacei;
1185 const label oldFacei = faceMap_[facei];
1186 map[patchFacei] = oldPatches[oldPatchi].whichFace(oldFacei);
1187 }
1188
1189 newBoundary.set
1190 (
1191 nNewPatches,
1192 oldPatches[oldPatchi].clone
1193 (
1194 subMeshPtr_().boundaryMesh(),
1195 nNewPatches,
1196 map,
1197 patchStart
1198 )
1199 );
1200 }
1201 else
1202 {
1203 // Patch still exists. Add it
1204 newBoundary.set
1205 (
1206 nNewPatches,
1207 oldPatches[oldPatchi].clone
1208 (
1209 subMeshPtr_().boundaryMesh(),
1210 nNewPatches,
1211 newSize,
1212 patchStart
1213 )
1214 );
1215 }
1216
1217 //Pout<< " " << oldPatches[oldPatchi].name() << " : "
1218 // << newSize << endl;
1219
1220 patchStart += newSize;
1221 patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1222 ++nNewPatches;
1223 }
1224
1225
1226 // Reset the patch lists
1227 newBoundary.resize(nNewPatches);
1228 patchMap_.resize(nNewPatches);
1229
1230
1231 // Add the fvPatches
1232 subMeshPtr_().addFvPatches(newBoundary, syncPar);
1233
1234 // Subset and add any zones
1235 subsetZones();
1236}
1237
1238
1240(
1241 const labelUList& selectedCells,
1242 const label patchID,
1243 const bool syncPar
1244)
1245{
1246 reset
1247 (
1248 BitSetOps::create(baseMesh().nCells(), selectedCells),
1249 patchID,
1250 syncPar
1251 );
1252}
1253
1254
1256(
1257 const labelHashSet& selectedCells,
1258 const label patchID,
1259 const bool syncPar
1260)
1261{
1262 reset
1263 (
1264 BitSetOps::create(baseMesh().nCells(), selectedCells),
1265 patchID,
1266 syncPar
1267 );
1268}
1269
1270
1272(
1273 const label regioni,
1274 const labelUList& regions,
1275 const label patchID,
1276 const bool syncPar
1277)
1278{
1279 reset
1280 (
1281 BitSetOps::create(baseMesh().nCells(), regioni, regions),
1282 patchID,
1283 syncPar
1284 );
1285}
1286
1287
1288// ************************************************************************* //
label n
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
Buffers for inter-processor communications streams (UOPstream, UIPstream).
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
void finishedSends(const bool wait=true)
Mark sends as done.
static void allGatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:103
A List obtained as a section of another List.
Definition: SubList.H:70
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:329
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
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
labelList sortedToc() const
The indices of the on bits as a sorted labelList.
Definition: bitSetI.H:533
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:521
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
Empty front and back plane patch. Used for 2-D geometries.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:80
bool checkHasSubMesh() const
FatalError if subset has not been performed.
Definition: fvMeshSubset.C:113
static word exposedPatchName
Name for exposed internal faces (default: oldInternalFaces)
Definition: fvMeshSubset.H:140
void clear()
Reset subMesh and all maps.
Definition: fvMeshSubset.C:481
void reset()
Reset subMesh and all maps. Same as clear()
Definition: fvMeshSubset.C:493
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
label nNonProcessor() const
The number of patches before the first processor patch.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
wordList names() const
Return a list of patch names.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
int myProcNo() const noexcept
Return processor number.
A class for handling words, derived from Foam::string.
Definition: word.H:68
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
patchWriters clear()
IOporosityModelList pZones(mesh)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
label nPoints
bitSet create(const label n, const labelHashSet &locations, const bool on=true)
Create a bitSet with length n with the specified on locations.
Definition: BitOps.C:212
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
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
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
static labelList subsetSubset(const label nElems, const labelUList &selectedElements, const labelUList &subsetMap)
Definition: fvMeshSubset.C:67
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:61
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
wordList patchNames(nPatches)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.