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