Loading...
Searching...
No Matches
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,2024 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#include "meshPointPatch.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
42 word fvMeshSubset::exposedPatchName("oldInternalFaces");
44}
45
46
47// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
48
49namespace
50{
51
52// Mark faces/points (with 0) in labelList
53inline void markUsed
54(
55 const Foam::labelUList& locations,
57)
58{
59 for (auto idx : locations)
60 {
61 map[idx] = 0;
62 }
63}
64
65} // End anonymous namespace
66
67
68namespace Foam
69{
70
71// Perform a subset of a subset
73(
74 const label nElems,
75 const labelUList& selectedElements, // First subset
76 const labelUList& subsetMap // Subset within first subset
77)
78{
79 if (selectedElements.empty() || subsetMap.empty())
80 {
81 // Trivial case
82 return labelList();
83 }
84
85 // Mark selected elements.
86 const bitSet selected(nElems, selectedElements);
87
88 // Count subset of selected elements
89 label n = 0;
90 forAll(subsetMap, i)
91 {
92 if (selected[subsetMap[i]])
93 {
94 ++n;
95 }
96 }
97
98 // Collect selected elements
99 labelList subsettedElements(n);
100 n = 0;
101
102 forAll(subsetMap, i)
103 {
104 if (selected[subsetMap[i]])
105 {
106 subsettedElements[n] = i;
107 ++n;
108 }
109 }
110
111 return subsettedElements;
113
114} // End namespace Foam
115
116
117// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
118
120{
121 if (!subMeshPtr_)
122 {
124 << "Mesh is not subsetted!" << nl
125 << abort(FatalError);
126
127 return false;
128 }
129
130 return true;
131}
132
133
134void Foam::fvMeshSubset::calcFaceFlipMap() const
135{
136 const labelList& subToBaseFace = faceMap();
137 const labelList& subToBaseCell = cellMap();
138
139 faceFlipMapPtr_.reset(new labelList(subToBaseFace.size()));
140 auto& faceFlipMap = *faceFlipMapPtr_;
141
142 // Only exposed internal faces might be flipped (since we don't do
143 // any cell renumbering, just compacting)
144 const label subInt = subMesh().nInternalFaces();
145
146 const labelList& subOwn = subMesh().faceOwner();
147 const labelList& own = baseMesh_.faceOwner();
148
149 for (label subFaceI = 0; subFaceI < subInt; ++subFaceI)
150 {
151 faceFlipMap[subFaceI] = subToBaseFace[subFaceI]+1;
152 }
153 for (label subFaceI = subInt; subFaceI < subOwn.size(); ++subFaceI)
154 {
155 const label faceI = subToBaseFace[subFaceI];
156 if (subToBaseCell[subOwn[subFaceI]] == own[faceI])
157 {
158 faceFlipMap[subFaceI] = faceI+1;
159 }
160 else
161 {
162 faceFlipMap[subFaceI] = -faceI-1;
163 }
164 }
165}
166
167
168void Foam::fvMeshSubset::doCoupledPatches
169(
170 const bool syncPar,
171 labelUList& nCellsUsingFace
172) const
173{
174 // Synchronize facesToSubset on both sides of coupled patches.
175 // Marks faces that become 'uncoupled' with 3.
176
177 const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
178
179 label nUncoupled = 0;
180
181 if (syncPar && UPstream::parRun())
182 {
184
185 // Send face usage across processor patches
186 if (!nCellsUsingFace.empty())
187 {
188 for (const polyPatch& pp : oldPatches)
189 {
190 const auto* procPatch = isA<processorPolyPatch>(pp);
191
192 if (procPatch)
193 {
194 const label nbrProci = procPatch->neighbProcNo();
196 UOPstream toNbr(nbrProci, pBufs);
197 toNbr <<
198 SubList<label>(nCellsUsingFace, pp.size(), pp.start());
199 }
200 }
201 }
202
203 pBufs.finishedSends();
204
205 // Receive face usage count and check for faces that become uncoupled.
206 for (const polyPatch& pp : oldPatches)
207 {
208 const auto* procPatch = isA<processorPolyPatch>(pp);
209
210 if (procPatch)
211 {
212 const label nbrProci = procPatch->neighbProcNo();
213
214 if (!pBufs.recvDataCount(nbrProci))
215 {
216 // Nothing to receive
217 continue;
218 }
219
220 labelList nbrCellsUsingFace;
222 UIPstream fromNbr(nbrProci, pBufs);
223 fromNbr >> nbrCellsUsingFace;
224 }
225
226 if (!nCellsUsingFace.empty() && !nbrCellsUsingFace.empty())
227 {
228 // Combine with this side.
229
230 forAll(pp, i)
231 {
232 if
233 (
234 nCellsUsingFace[pp.start()+i] == 1
235 && nbrCellsUsingFace[i] == 0
236 )
237 {
238 // Face's neighbour is no longer there.
239 // Mark face off as coupled
240 nCellsUsingFace[pp.start()+i] = 3;
241 ++nUncoupled;
242 }
243 }
244 }
245 }
246 }
247 }
248
249 // Do same for cyclics.
250 for (const polyPatch& pp : oldPatches)
251 {
253
254 if (cpp && !nCellsUsingFace.empty())
255 {
256 const auto& cycPatch = *cpp;
257
258 forAll(cycPatch, i)
259 {
260 label thisFacei = cycPatch.start() + i;
261 label otherFacei = cycPatch.transformGlobalFace(thisFacei);
262
263 if
264 (
265 nCellsUsingFace[thisFacei] == 1
266 && nCellsUsingFace[otherFacei] == 0
267 )
268 {
269 nCellsUsingFace[thisFacei] = 3;
270 ++nUncoupled;
271 }
272 }
273 }
274 }
275
276 if (syncPar)
277 {
278 reduce(nUncoupled, sumOp<label>());
279 }
280
281 if (nUncoupled)
282 {
283 Info<< "Uncoupled " << nUncoupled << " faces on coupled patches. "
284 << "(processorPolyPatch, cyclicPolyPatch)" << endl;
285 }
286}
287
288
289void Foam::fvMeshSubset::subsetZones()
290{
291 // Keep all zones, even if zero size.
292
293 #ifdef FULLDEBUG
294 checkHasSubMesh();
295 #endif
296
297 auto& newSubMesh = subMeshPtr_();
298
299 // PointZones
300
301 const pointZoneMesh& pointZones = baseMesh().pointZones();
302
303 List<pointZone*> pZones(pointZones.size());
304
305 forAll(pointZones, zonei)
306 {
307 const pointZone& pz = pointZones[zonei];
308
309 pZones[zonei] = new pointZone
310 (
311 pz.name(),
312 subsetSubset(baseMesh().nPoints(), pz, pointMap()),
313 zonei,
314 newSubMesh.pointZones()
315 );
316 }
317
318
319 // FaceZones
320 // Do we need to remove zones where the side we're interested in
321 // no longer exists? Guess not.
322
323 const faceZoneMesh& faceZones = baseMesh().faceZones();
324
325 List<faceZone*> fZones(faceZones.size());
326
327 forAll(faceZones, zonei)
328 {
329 const faceZone& fz = faceZones[zonei];
330
331 // Expand faceZone to full mesh
332 // +1 : part of faceZone, flipped
333 // -1 : ,, , unflipped
334 // 0 : not part of faceZone
335 labelList zone(baseMesh().nFaces(), Zero);
336 forAll(fz, j)
337 {
338 zone[fz[j]] = (fz.flipMap()[j] ? 1 : -1);
339 }
340
341 // Select faces
342 label nSub = 0;
343 forAll(faceMap(), j)
344 {
345 if (zone[faceMap()[j]] != 0)
346 {
347 ++nSub;
348 }
349 }
350 labelList subAddressing(nSub);
351 boolList subFlipStatus(nSub);
352 nSub = 0;
353 forAll(faceMap(), subFacei)
354 {
355 const label meshFacei = faceMap()[subFacei];
356 if (zone[meshFacei] != 0)
357 {
358 subAddressing[nSub] = subFacei;
359 const label subOwner = subMesh().faceOwner()[subFacei];
360 const label baseOwner = baseMesh().faceOwner()[meshFacei];
361 // If subowner is the same cell as the base keep the flip status
362 const bool sameOwner = (cellMap()[subOwner] == baseOwner);
363 const bool flip = (zone[meshFacei] == 1);
364 subFlipStatus[nSub] = (sameOwner == flip);
365
366 ++nSub;
367 }
368 }
369
370 fZones[zonei] = new faceZone
371 (
372 fz.name(),
373 subAddressing,
374 subFlipStatus,
375 zonei,
376 newSubMesh.faceZones()
377 );
378 }
379
380
381 // Cell Zones
382
383 const cellZoneMesh& cellZones = baseMesh().cellZones();
384
385 List<cellZone*> cZones(cellZones.size());
386
387 forAll(cellZones, zonei)
388 {
389 const cellZone& cz = cellZones[zonei];
390
391 cZones[zonei] = new cellZone
392 (
393 cz.name(),
394 subsetSubset(baseMesh().nCells(), cz, cellMap()),
395 zonei,
396 newSubMesh.cellZones()
397 );
398 }
399
400 // Add the zones
401 newSubMesh.addZones(pZones, fZones, cZones);
402}
403
404
405// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
408:
409 baseMesh_(baseMesh)
410{}
411
412
416{
417 reset(Foam::zero{});
418}
419
420
422(
423 const fvMesh& baseMesh,
424 const bitSet& selectedCells,
425 const label patchID,
426 const bool syncPar
427)
430{
431 reset(selectedCells, patchID, syncPar);
432}
433
434
436(
437 const fvMesh& baseMesh,
438 const labelUList& selectedCells,
439 const label patchID,
440 const bool syncPar
441)
444{
445 reset(selectedCells, patchID, syncPar);
446}
447
448
450(
451 const fvMesh& baseMesh,
452 const labelHashSet& selectedCells,
453 const label patchID,
454 const bool syncPar
455)
458{
459 reset(selectedCells, patchID, syncPar);
460}
461
462
464(
465 const fvMesh& baseMesh,
466 const label regioni,
467 const labelUList& regions,
468 const label patchID,
469 const bool syncPar
470)
471:
472 fvMeshSubset(baseMesh)
474 reset(regioni, regions, patchID, syncPar);
475}
476
477
478// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
479
481{
482 subMeshPtr_.reset(nullptr);
483 faceFlipMapPtr_.reset(nullptr);
484
485 pointMap_.clear();
486 faceMap_.clear();
487 cellMap_.clear();
488 patchMap_.clear();
489}
490
493{
494 clear();
495}
496
497
499(
500 autoPtr<fvMesh>&& subMeshPtr,
501 labelList&& pointMap,
503 labelList&& cellMap,
504 labelList&& patchMap
505)
506{
507 subMeshPtr_.reset(std::move(subMeshPtr));
508 faceFlipMapPtr_.reset(nullptr);
509
510 pointMap_ = std::move(pointMap);
511 faceMap_ = std::move(faceMap);
512 cellMap_ = std::move(cellMap);
513 patchMap_ = std::move(patchMap);
514
515 // Sanity
516 if (!subMeshPtr_)
517 {
518 clear();
519 }
520}
521
522
524{
525 // Was old pointMesh present?
526 const auto* basePointMeshPtr =
527 baseMesh_.thisDb().cfindObject<pointMesh>(pointMesh::typeName);
528 if (basePointMeshPtr)
529 {
530 DebugPout<< "fvMeshSubset::reset(Foam::zero) :"
531 << " Detected pointMesh" << endl;
532 }
533
534
535 clear();
536
537 // Create zero-sized subMesh
538 subMeshPtr_.reset
539 (
540 new fvMesh
541 (
543 (
544 baseMesh_.name(),
545 baseMesh_.time().timeName(),
546 baseMesh_.time(),
547 IOobject::NO_READ, // Do not read any dictionaries
549 ),
550 baseMesh_, // Get dictionaries from base mesh
551 Foam::zero{} // zero-sized
552 // Uses syncPar (bounds) - should generally be OK
553 )
554 );
555 auto& newSubMesh = subMeshPtr_();
556
557
558 // Clone non-processor patches
559 {
560 const polyBoundaryMesh& oldBoundary = baseMesh_.boundaryMesh();
561 const polyBoundaryMesh& newBoundary = newSubMesh.boundaryMesh();
562
563 polyPatchList newPatches(oldBoundary.nNonProcessor());
564
565 patchMap_ = identity(newPatches.size());
566
567 forAll(newPatches, patchi)
568 {
569 newPatches.set
570 (
571 patchi,
572 oldBoundary[patchi].clone
573 (
574 newBoundary,
575 patchi,
576 0, // patch size
577 0 // patch start
578 )
579 );
580 }
581
582 // Add patches - make sure we don't trigger any parallel side effects
583 newSubMesh.addFvPatches(newPatches, false);
584 }
585
586
587 // Clone old additional point patches
588 if (basePointMeshPtr)
589 {
590 DebugPout<< "Subsetting pointMesh" << endl;
591 const auto& basePointMesh = *basePointMeshPtr;
592 const auto& oldPointBoundary = basePointMesh.boundary();
593
594 // 1. Generate pointBoundaryMesh from polyBoundaryMesh (so ignoring
595 // any additional patches
596 const auto& newSubPointMesh = pointMesh::New(newSubMesh);
597
598 auto& newBoundary =
599 const_cast<pointBoundaryMesh&>(newSubPointMesh.boundary());
600
601 // Start off from (poly)patch map
602 pointPatchMap_ = patchMap_;
603
604 // 2. Explicitly add subsetted meshPointPatches
605 for (const auto& oldPointPatch : oldPointBoundary)
606 {
607 const auto* mppPtr = isA<meshPointPatch>(oldPointPatch);
608 if (mppPtr && (newBoundary.findPatchID(mppPtr->name()) == -1))
609 {
610 newBoundary.push_back
611 (
612 mppPtr->clone
613 (
614 newBoundary,
615 newBoundary.size(),
616 labelList::null(), // map
617 labelList::null() // map
618 )
619 );
620 }
621 }
622
623 // Extend patchMap with -1
624 pointPatchMap_.setSize(newBoundary.size(), -1);
626
627 // Add the zones
628 subsetZones();
629}
630
631
633(
634 const bitSet& selectedCells,
635 const label patchID,
636 const bool syncPar
637)
638{
639 // Was old pointMesh present?
640 const auto* basePointMeshPtr =
641 baseMesh_.thisDb().cfindObject<pointMesh>(pointMesh::typeName);
642 if (basePointMeshPtr)
643 {
644 DebugPout<< "fvMeshSubset::reset(const bitSet&) :"
645 << " Detected pointMesh" << endl;
646 }
647
648
649 // Clear all old maps and pointers
650 clear();
651
652 const cellList& oldCells = baseMesh().cells();
653 const faceList& oldFaces = baseMesh().faces();
654 const pointField& oldPoints = baseMesh().points();
655 const labelList& oldOwner = baseMesh().faceOwner();
656 const labelList& oldNeighbour = baseMesh().faceNeighbour();
657 const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
658 const label oldNInternalFaces = baseMesh().nInternalFaces();
659
660 // Initial check on patches before doing anything time consuming.
661
662 label wantedPatchID = patchID;
663
664 if (wantedPatchID == -1)
665 {
666 // No explicit patch specified. Put in oldInternalFaces patch.
667 // Check if patch with this name already exists.
668 wantedPatchID = oldPatches.findPatchID(exposedPatchName);
669 }
670 else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
671 {
673 << "Non-existing patch index " << wantedPatchID << endl
674 << "Should be between 0 and " << oldPatches.size()-1
675 << abort(FatalError);
676 }
677
678 // The selected cells - sorted in ascending order
679 cellMap_ = selectedCells.sortedToc();
680
681 // The selectedCells should normally be in the [0,nCells) range,
682 // but it is better not to trust that.
683 {
684 label len = cellMap_.size();
685 for
686 (
687 label i = len-1;
688 i >= 0 && (cellMap_[i] >= oldCells.size());
689 --i
690 )
691 {
692 len = i;
693 }
694 cellMap_.resize(len);
695 }
696
697
698 // Mark all used faces. Count number of cells using them
699 // 0: face not used anymore
700 // 1: face used by one cell, face becomes/stays boundary face
701 // 2: face still used and remains internal face
702 // 3: face coupled and used by one cell only (so should become normal,
703 // non-coupled patch face)
704 //
705 // Note that this is not really necessary - but means we can size things
706 // correctly. Also makes handling coupled faces much easier.
707
708 labelList nCellsUsingFace(oldFaces.size(), Zero);
709
710 label nFacesInSet = 0;
711
712 forAll(oldFaces, oldFacei)
713 {
714 bool faceUsed = false;
715
716 if (selectedCells.test(oldOwner[oldFacei]))
717 {
718 ++nCellsUsingFace[oldFacei];
719 faceUsed = true;
720 }
721
722 if
723 (
724 baseMesh().isInternalFace(oldFacei)
725 && selectedCells.test(oldNeighbour[oldFacei])
726 )
727 {
728 ++nCellsUsingFace[oldFacei];
729 faceUsed = true;
730 }
731
732 if (faceUsed)
733 {
734 ++nFacesInSet;
735 }
736 }
737 faceMap_.resize(nFacesInSet);
738
739 // Handle coupled faces. Modifies patch faces to be uncoupled to 3.
740 doCoupledPatches(syncPar, nCellsUsingFace);
741
742
743 // See which patch to use for exposed internal faces.
744 label oldInternalPatchID = 0;
745
746 // Insert faces before which patch
747 label nextPatchID = oldPatches.size();
748
749 // old to new patches
750 labelList globalPatchMap(oldPatches.size());
751
752 // New patch size
753 label nbSize = oldPatches.size();
754
755 if (wantedPatchID == -1)
756 {
757 // Create 'oldInternalFaces' patch at the end (or before
758 // processorPatches)
759 // and put all exposed internal faces in there.
760
761 forAll(oldPatches, patchi)
762 {
763 if (isA<processorPolyPatch>(oldPatches[patchi]))
764 {
765 nextPatchID = patchi;
766 break;
767 }
768 ++oldInternalPatchID;
769 }
770
771 ++nbSize;
772
773 // adapt old to new patches for inserted patch
774 for (label oldPatchi = 0; oldPatchi < nextPatchID; oldPatchi++)
775 {
776 globalPatchMap[oldPatchi] = oldPatchi;
777 }
778 for
779 (
780 label oldPatchi = nextPatchID;
781 oldPatchi < oldPatches.size();
782 oldPatchi++
783 )
784 {
785 globalPatchMap[oldPatchi] = oldPatchi+1;
786 }
787 }
788 else
789 {
790 oldInternalPatchID = wantedPatchID;
791 nextPatchID = wantedPatchID+1;
792
793 // old to new patches
794 globalPatchMap = identity(oldPatches.size());
795 }
796
797 labelList boundaryPatchSizes(nbSize, Zero);
798
799
800 // Make a global-to-local point map
801 labelList globalPointMap(oldPoints.size(), -1);
802 labelList globalFaceMap(oldFaces.size(), -1);
803
804 label facei = 0;
805
806 // 1. Pick up all preserved internal faces.
807 for (label oldFacei = 0; oldFacei < oldNInternalFaces; ++oldFacei)
808 {
809 if (nCellsUsingFace[oldFacei] == 2)
810 {
811 globalFaceMap[oldFacei] = facei;
812 faceMap_[facei++] = oldFacei;
813
814 // Mark all points from the face
815 markUsed(oldFaces[oldFacei], globalPointMap);
816 }
817 }
818
819 // These are all the internal faces in the mesh.
820 const label nInternalFaces = facei;
821
822 // 2. Boundary faces up to where we want to insert old internal faces
823 for
824 (
825 label oldPatchi = 0;
826 oldPatchi < oldPatches.size()
827 && oldPatchi < nextPatchID;
828 oldPatchi++
829 )
830 {
831 const polyPatch& oldPatch = oldPatches[oldPatchi];
832
833 label oldFacei = oldPatch.start();
834
835 forAll(oldPatch, i)
836 {
837 if (nCellsUsingFace[oldFacei] == 1)
838 {
839 // Boundary face is kept.
840
841 // Mark face and increment number of points in set
842 globalFaceMap[oldFacei] = facei;
843 faceMap_[facei++] = oldFacei;
844
845 // Mark all points from the face
846 markUsed(oldFaces[oldFacei], globalPointMap);
847
848 // Increment number of patch faces
849 ++boundaryPatchSizes[globalPatchMap[oldPatchi]];
850 }
851 ++oldFacei;
852 }
853 }
854
855 // 3a. old internal faces that have become exposed.
856 for (label oldFacei = 0; oldFacei < oldNInternalFaces; ++oldFacei)
857 {
858 if (nCellsUsingFace[oldFacei] == 1)
859 {
860 globalFaceMap[oldFacei] = facei;
861 faceMap_[facei++] = oldFacei;
862
863 // Mark all points from the face
864 markUsed(oldFaces[oldFacei], globalPointMap);
865
866 // Increment number of patch faces
867 ++boundaryPatchSizes[oldInternalPatchID];
868 }
869 }
870
871 // 3b. coupled patch faces that have become uncoupled.
872 for
873 (
874 label oldFacei = oldNInternalFaces;
875 oldFacei < oldFaces.size();
876 oldFacei++
877 )
878 {
879 if (nCellsUsingFace[oldFacei] == 3)
880 {
881 globalFaceMap[oldFacei] = facei;
882 faceMap_[facei++] = oldFacei;
883
884 // Mark all points from the face
885 markUsed(oldFaces[oldFacei], globalPointMap);
886
887 // Increment number of patch faces
888 ++boundaryPatchSizes[oldInternalPatchID];
889 }
890 }
891
892 // 4. Remaining boundary faces
893 for
894 (
895 label oldPatchi = nextPatchID;
896 oldPatchi < oldPatches.size();
897 oldPatchi++
898 )
899 {
900 const polyPatch& oldPatch = oldPatches[oldPatchi];
901
902 label oldFacei = oldPatch.start();
903
904 forAll(oldPatch, i)
905 {
906 if (nCellsUsingFace[oldFacei] == 1)
907 {
908 // Boundary face is kept.
909
910 // Mark face and increment number of points in set
911 globalFaceMap[oldFacei] = facei;
912 faceMap_[facei++] = oldFacei;
913
914 // Mark all points from the face
915 markUsed(oldFaces[oldFacei], globalPointMap);
916
917 // Increment number of patch faces
918 ++boundaryPatchSizes[globalPatchMap[oldPatchi]];
919 }
920 ++oldFacei;
921 }
922 }
923
924 if (facei != nFacesInSet)
925 {
927 << "Problem" << abort(FatalError);
928 }
929
930
931 // Grab the points map
932 label nPointsInSet = 0;
933
934 forAll(globalPointMap, pointi)
935 {
936 if (globalPointMap[pointi] != -1)
937 {
938 ++nPointsInSet;
939 }
940 }
941 pointMap_.setSize(nPointsInSet);
942
943 nPointsInSet = 0;
944
945 forAll(globalPointMap, pointi)
946 {
947 if (globalPointMap[pointi] != -1)
948 {
949 pointMap_[nPointsInSet] = pointi;
950 globalPointMap[pointi] = nPointsInSet;
951 ++nPointsInSet;
952 }
953 }
954
955
956 //Pout<< "Number of points,faces,cells in new mesh : "
957 // << pointMap_.size() << ' '
958 // << faceMap_.size() << ' '
959 // << cellMap_.size() << nl;
960
961
962 // Make a new mesh
963
964 //
965 // Create new points
966 //
967 pointField newPoints(oldPoints, pointMap_);
968
969
970 //
971 // Create new faces
972 //
973
974 faceList newFaces(faceMap_.size());
975 {
976 auto iter = newFaces.begin();
977 const auto& renumbering = globalPointMap;
978
979 // Internal faces
980 for (label facei = 0; facei < nInternalFaces; ++facei)
981 {
982 face& newItem = *iter;
983 ++iter;
984
985 const face& oldItem = oldFaces[faceMap_[facei]];
986
987 // Copy relabelled
988 newItem.resize(oldItem.size());
989
990 forAll(oldItem, i)
991 {
992 newItem[i] = renumbering[oldItem[i]];
993 }
994 }
995
996 // Boundary faces - may need to be flipped
997 for (label facei = nInternalFaces; facei < faceMap_.size(); ++facei)
998 {
999 const label oldFacei = faceMap_[facei];
1000
1001 face& newItem = *iter;
1002 ++iter;
1003
1004 const face oldItem =
1005 (
1006 (
1007 baseMesh().isInternalFace(oldFacei)
1008 && selectedCells.test(oldNeighbour[oldFacei])
1009 && !selectedCells.test(oldOwner[oldFacei])
1010 )
1011 // Face flipped to point outwards
1012 ? oldFaces[oldFacei].reverseFace()
1013 : oldFaces[oldFacei]
1014 );
1015
1016 // Copy relabelled
1017 newItem.resize(oldItem.size());
1018
1019 forAll(oldItem, i)
1020 {
1021 newItem[i] = renumbering[oldItem[i]];
1022 }
1023 }
1024 }
1025
1026
1027 //
1028 // Create new cells
1029 //
1030
1031 cellList newCells(cellMap_.size());
1032 {
1033 auto iter = newCells.begin();
1034 const auto& renumbering = globalFaceMap;
1035
1036 for (const label oldCelli : cellMap_)
1037 {
1038 cell& newItem = *iter;
1039 ++iter;
1040
1041 const labelList& oldItem = oldCells[oldCelli];
1042
1043 // Copy relabelled
1044 newItem.resize(oldItem.size());
1045
1046 forAll(oldItem, i)
1047 {
1048 newItem[i] = renumbering[oldItem[i]];
1049 }
1050 }
1051 }
1052
1053
1054 // Make a new mesh
1055 //
1056 // Note that mesh gets registered with same name as original mesh.
1057 // This is not proper but cannot be avoided since otherwise
1058 // surfaceInterpolation cannot find its fvSchemes.
1059 // It will try to read, for example. "system/region0SubSet/fvSchemes"
1060 //
1061 subMeshPtr_ = autoPtr<fvMesh>::New
1062 (
1063 IOobject
1064 (
1065 baseMesh_.name(),
1066 baseMesh_.time().timeName(),
1067 baseMesh_.time(),
1068 IOobject::NO_READ, // Do not read any dictionaries
1070 ),
1071 baseMesh_, // Get dictionaries from base mesh
1072 std::move(newPoints),
1073 std::move(newFaces),
1074 std::move(newCells),
1075 syncPar // parallel synchronisation
1076 );
1077
1078
1079 // Add old patches
1080 polyPatchList newBoundary(nbSize);
1081 patchMap_.resize(nbSize);
1082 label nNewPatches = 0;
1083 label patchStart = nInternalFaces;
1084
1085
1086 // For parallel: only remove patch if none of the processors has it.
1087 // This only gets done for patches before the one being inserted
1088 // (so patches < nextPatchID)
1089
1090 // Get sum of patch sizes. Zero if patch can be deleted.
1091 labelList globalPatchSizes(boundaryPatchSizes);
1092 globalPatchSizes.setSize(nextPatchID);
1093
1094 if (syncPar && Pstream::parRun())
1095 {
1096 // Get patch names (up to nextPatchID)
1098 patchNames[Pstream::myProcNo()] = oldPatches.names();
1099 patchNames[Pstream::myProcNo()].setSize(nextPatchID);
1101
1102 // Get patch sizes (up to nextPatchID).
1103 // Note that up to nextPatchID the globalPatchMap is an identity so
1104 // no need to index through that.
1105 Pstream::listReduce(globalPatchSizes, sumOp<label>());
1106
1107 // Now all processors have all the patchnames.
1108 // Decide: if all processors have the same patch names and size is zero
1109 // everywhere remove the patch.
1110 bool samePatches = true;
1111
1112 for (label proci = 1; proci < patchNames.size(); ++proci)
1113 {
1114 if (patchNames[proci] != patchNames[0])
1115 {
1116 samePatches = false;
1117 break;
1118 }
1119 }
1120
1121 if (!samePatches)
1122 {
1123 // Patchnames not sync on all processors so disable removal of
1124 // zero sized patches.
1125 globalPatchSizes = labelMax;
1126 }
1127 }
1128
1129
1130 // Old patches
1131
1132 for
1133 (
1134 label oldPatchi = 0;
1135 oldPatchi < oldPatches.size()
1136 && oldPatchi < nextPatchID;
1137 oldPatchi++
1138 )
1139 {
1140 const label newSize = boundaryPatchSizes[globalPatchMap[oldPatchi]];
1141
1142 if (oldInternalPatchID != oldPatchi)
1143 {
1144 // Pure subset of patch faces (no internal faces added to this
1145 // patch). Can use mapping.
1146 labelList map(newSize);
1147 for (label patchFacei = 0; patchFacei < newSize; patchFacei++)
1148 {
1149 const label facei = patchStart+patchFacei;
1150 const label oldFacei = faceMap_[facei];
1151 map[patchFacei] = oldPatches[oldPatchi].whichFace(oldFacei);
1152 }
1153
1154 newBoundary.set
1155 (
1156 nNewPatches,
1157 oldPatches[oldPatchi].clone
1158 (
1159 subMeshPtr_().boundaryMesh(),
1160 nNewPatches,
1161 map,
1162 patchStart
1163 )
1164 );
1165 }
1166 else
1167 {
1168 // Clone (even if 0 size)
1169 newBoundary.set
1170 (
1171 nNewPatches,
1172 oldPatches[oldPatchi].clone
1173 (
1174 subMeshPtr_().boundaryMesh(),
1175 nNewPatches,
1176 newSize,
1177 patchStart
1178 )
1179 );
1180 }
1181
1182 patchStart += newSize;
1183 patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1184 ++nNewPatches;
1185 }
1186
1187 // Inserted patch
1188
1189 label newInternalPatchID = -1;
1190
1191 if (wantedPatchID == -1)
1192 {
1193 label oldInternalSize = boundaryPatchSizes[oldInternalPatchID];
1194
1195 if (syncPar)
1196 {
1197 reduce(oldInternalSize, sumOp<label>());
1198 }
1199
1200 // Newly created patch so is at end. Check if any faces in it.
1201 if (oldInternalSize > 0)
1202 {
1203 newBoundary.set
1204 (
1205 nNewPatches,
1206 new emptyPolyPatch
1207 (
1208 exposedPatchName,
1209 boundaryPatchSizes[oldInternalPatchID],
1210 patchStart,
1211 nNewPatches,
1212 subMeshPtr_().boundaryMesh(),
1213 emptyPolyPatch::typeName
1214 )
1215 );
1216
1217 //Pout<< " " << exposedPatchName << " : "
1218 // << boundaryPatchSizes[oldInternalPatchID] << endl;
1219
1220 // The index for the first patch is -1 as it originates from
1221 // the internal faces
1222 patchStart += boundaryPatchSizes[oldInternalPatchID];
1223 patchMap_[nNewPatches] = -1;
1224 newInternalPatchID = nNewPatches;
1225 ++nNewPatches;
1226 }
1227 }
1228
1229 // Old patches
1230
1231 for
1232 (
1233 label oldPatchi = nextPatchID;
1234 oldPatchi < oldPatches.size();
1235 oldPatchi++
1236 )
1237 {
1238 const label newSize = boundaryPatchSizes[globalPatchMap[oldPatchi]];
1239
1240 if (oldInternalPatchID != oldPatchi)
1241 {
1242 // Pure subset of patch faces (no internal faces added to this
1243 // patch). Can use mapping.
1244 labelList map(newSize);
1245 for (label patchFacei = 0; patchFacei < newSize; patchFacei++)
1246 {
1247 const label facei = patchStart+patchFacei;
1248 const label oldFacei = faceMap_[facei];
1249 map[patchFacei] = oldPatches[oldPatchi].whichFace(oldFacei);
1250 }
1251
1252 newBoundary.set
1253 (
1254 nNewPatches,
1255 oldPatches[oldPatchi].clone
1256 (
1257 subMeshPtr_().boundaryMesh(),
1258 nNewPatches,
1259 map,
1260 patchStart
1261 )
1262 );
1263 }
1264 else
1265 {
1266 // Patch still exists. Add it
1267 newBoundary.set
1268 (
1269 nNewPatches,
1270 oldPatches[oldPatchi].clone
1271 (
1272 subMeshPtr_().boundaryMesh(),
1273 nNewPatches,
1274 newSize,
1275 patchStart
1276 )
1277 );
1278 }
1279
1280 //Pout<< " " << oldPatches[oldPatchi].name() << " : "
1281 // << newSize << endl;
1282
1283 patchStart += newSize;
1284 patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1285 ++nNewPatches;
1286 }
1287
1288
1289 // Reset the patch lists
1290 newBoundary.resize(nNewPatches);
1291 patchMap_.resize(nNewPatches);
1292
1293
1294 // Add the fvPatches
1295 subMeshPtr_().addFvPatches(newBoundary, syncPar);
1296
1297 // Subset and add any zones
1298 subsetZones();
1299
1300
1301 if (basePointMeshPtr)
1302 {
1303 DebugPout<< "Subsetting pointMesh" << endl;
1304 const auto& basePointMesh = *basePointMeshPtr;
1305 const auto& oldPointBoundary = basePointMesh.boundary();
1306
1307 // 1. Generate pointBoundaryMesh from polyBoundaryMesh (so ignoring
1308 // any additional patches
1309 const auto& newSubPointMesh = pointMesh::New(subMeshPtr_());
1310
1311 pointPatchMap_ = patchMap_;
1312
1313 auto& newBoundary =
1314 const_cast<pointBoundaryMesh&>(newSubPointMesh.boundary());
1315
1316
1317 // 2. Explicitly add subsetted meshPointPatches
1318 labelList oldToNewPoints(baseMesh_.nPoints(), -1);
1319 forAll(pointMap_, i)
1320 {
1321 oldToNewPoints[pointMap_[i]] = i;
1322 }
1323
1324
1325 // Add meshPointPatches
1326 pointPatchMap_.setSize(newBoundary.size(), -1);
1327
1328 for (const auto& oldPointPatch : oldPointBoundary)
1329 {
1330 const auto* mppPtr = isA<meshPointPatch>(oldPointPatch);
1331 if (mppPtr && (newBoundary.findPatchID(mppPtr->name()) == -1))
1332 {
1333 const auto& mp = mppPtr->meshPoints();
1334 DynamicList<label> subPointMap(mp.size());
1335 forAll(mp, i)
1336 {
1337 const label newPointi = oldToNewPoints[mp[i]];
1338 if (newPointi != -1)
1339 {
1340 subPointMap.append(i);
1341 }
1342 }
1343
1344 pointPatchMap_.push_back(mppPtr->index());
1345
1346 newBoundary.push_back
1347 (
1348 mppPtr->clone
1349 (
1350 newBoundary,
1351 newBoundary.size(),
1352 subPointMap, // map
1353 oldToNewPoints
1354 )
1355 );
1356 }
1357 }
1358
1359
1360 // 3. rotate into place:
1361 // - global patches (including meshPointPatches)
1362 // - optional 'internalFaces' patch
1363 // - processor patches
1364 labelList oldToNew(newBoundary.size());
1365 label newPatchi = 0;
1366 forAll(newBoundary, patchi)
1367 {
1368 if
1369 (
1370 patchi != newInternalPatchID
1371 && !isA<processorPointPatch>(newBoundary[patchi])
1372 )
1373 {
1374 oldToNew[patchi] = newPatchi++;
1375 }
1376 }
1377 if (newInternalPatchID != -1)
1378 {
1379 oldToNew[newInternalPatchID] = newPatchi++;
1380 }
1381 forAll(newBoundary, patchi)
1382 {
1383 if (isA<processorPointPatch>(newBoundary[patchi]))
1384 {
1385 oldToNew[patchi] = newPatchi++;
1386 }
1388 newBoundary.reorder(oldToNew, true);
1389 inplaceReorder(oldToNew, pointPatchMap_);
1390 }
1391}
1392
1393
1395(
1396 const labelUList& selectedCells,
1397 const label patchID,
1398 const bool syncPar
1399)
1400{
1401 reset
1402 (
1403 BitSetOps::create(baseMesh().nCells(), selectedCells),
1404 patchID,
1405 syncPar
1406 );
1407}
1408
1409
1411(
1412 const labelHashSet& selectedCells,
1413 const label patchID,
1414 const bool syncPar
1415)
1416{
1417 reset
1418 (
1419 BitSetOps::create(baseMesh().nCells(), selectedCells),
1420 patchID,
1421 syncPar
1422 );
1423}
1424
1425
1427(
1428 const label regioni,
1429 const labelUList& regions,
1430 const label patchID,
1431 const bool syncPar
1432)
1433{
1434 reset
1435 (
1436 BitSetOps::create(baseMesh().nCells(), regioni, regions),
1437 patchID,
1438 syncPar
1439 );
1440}
1441
1442
1443// ************************************************************************* //
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
static const List< label > & null() noexcept
Definition List.H:138
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
label recvDataCount(const label proci) const
Number of unconsumed receive bytes for the specified processor. Must call finishedSends() or other fi...
void finishedSends(const bool wait=true)
Mark the send phase as being finished.
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void setSize(const label n)
Same as resize().
Definition PtrList.H:357
void push_back(T *ptr)
Append an element to the end of the list.
Definition PtrListI.H:131
void resize(const label newLen)
Adjust size of PtrList.
Definition PtrList.C:124
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UIPstream.H:313
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UOPstream.H:408
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
void reorder(const labelUList &oldToNew, const bool check=false)
Reorder elements. Reordering must be unique (ie, shuffle).
Definition UPtrList.C:62
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
labelList sortedToc() const
The indices of the on bits as a sorted labelList.
Definition bitSetI.H:441
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition bitSet.H:334
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A subset of mesh cells.
Definition cellZone.H:61
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
Cyclic plane patch.
Empty front and back plane patch. Used for 2-D geometries.
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
fvMeshSubset(const fvMeshSubset &)=delete
No copy construct.
static word exposedPatchName
Name for exposed internal faces (default: oldInternalFaces).
const fvMesh & baseMesh() const noexcept
Original mesh.
const labelList & faceMap() const
Return face map.
const labelList & cellMap() const
Return cell map.
bool checkHasSubMesh() const
FatalError if subset has not been performed.
const labelList & patchMap() const
Return patch map.
const labelList & pointMap() const
Return point map.
void clear()
Reset subMesh and all maps.
void reset()
Reset subMesh and all maps. Same as clear().
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A pointBoundaryMesh is a pointPatch list with registered IO, a reference to the associated pointMesh,...
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
A subset of mesh points.
Definition pointZone.H:64
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
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:73
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
Base class for mesh zones.
Definition zone.H:63
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
limits reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL))
IOporosityModelList pZones(mesh)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
label nPoints
surface1 clear()
#define DebugPout
Report an information message using Foam::Pout.
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:233
const dimensionedScalar mp
Proton mass.
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with pointZone content on a polyMesh.
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
static labelList subsetSubset(const label nElems, const labelUList &selectedElements, const labelUList &subsetMap)
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
constexpr label labelMax
Definition label.H:55
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const direction noexcept
Definition scalarImpl.H:265
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
label newPointi
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299