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-2020 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_)
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_.reset(new labelList(subToBaseFace.size()));
91  auto& faceFlipMap = *faceFlipMapPtr_;
92 
93  // Only exposed internal faces might be flipped (since we don't do
94  // any cell renumbering, just compacting)
95  const label subInt = subMesh().nInternalFaces();
96 
97  const labelList& subOwn = subMesh().faceOwner();
98  const labelList& own = baseMesh_.faceOwner();
99 
100  for (label subFaceI = 0; subFaceI < subInt; ++subFaceI)
101  {
102  faceFlipMap[subFaceI] = subToBaseFace[subFaceI]+1;
103  }
104  for (label subFaceI = subInt; subFaceI < subOwn.size(); ++subFaceI)
105  {
106  const label faceI = subToBaseFace[subFaceI];
107  if (subToBaseCell[subOwn[subFaceI]] == own[faceI])
108  {
109  faceFlipMap[subFaceI] = faceI+1;
110  }
111  else
112  {
113  faceFlipMap[subFaceI] = -faceI-1;
114  }
115  }
116 }
117 
118 
119 void Foam::fvMeshSubset::doCoupledPatches
120 (
121  const bool syncPar,
122  labelList& nCellsUsingFace
123 ) const
124 {
125  // Synchronize facesToSubset on both sides of coupled patches.
126  // Marks faces that become 'uncoupled' with 3.
127 
128  const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
129 
130  label nUncoupled = 0;
131 
132  if (syncPar && Pstream::parRun())
133  {
134  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
135 
136  // Send face usage across processor patches
137  for (const polyPatch& pp : oldPatches)
138  {
139  if (isA<processorPolyPatch>(pp))
140  {
141  const processorPolyPatch& procPatch =
142  refCast<const processorPolyPatch>(pp);
143 
144  UOPstream toNeighbour(procPatch.neighbProcNo(), pBufs);
145 
146  if (!nCellsUsingFace.empty())
147  {
148  toNeighbour <<
149  SubList<label>(nCellsUsingFace, pp.size(), pp.start());
150  }
151  else
152  {
153  toNeighbour << labelList();
154  }
155  }
156  }
157 
158  pBufs.finishedSends();
159 
160  // Receive face usage count and check for faces that become uncoupled.
161  for (const polyPatch& pp : oldPatches)
162  {
163  if (isA<processorPolyPatch>(pp))
164  {
165  const processorPolyPatch& procPatch =
166  refCast<const processorPolyPatch>(pp);
167 
168  UIPstream fromNeighbour(procPatch.neighbProcNo(), pBufs);
169 
170  const labelList nbrList(fromNeighbour);
171 
172  // Combine with this side.
173 
174  if (!nCellsUsingFace.empty())
175  {
176  const labelList& nbrCellsUsingFace(nbrList);
177 
178  // Combine with this side.
179 
180  forAll(pp, i)
181  {
182  if
183  (
184  nCellsUsingFace[pp.start()+i] == 1
185  && nbrCellsUsingFace[i] == 0
186  )
187  {
188  // Face's neighbour is no longer there. Mark face
189  // off as coupled
190  nCellsUsingFace[pp.start()+i] = 3;
191  ++nUncoupled;
192  }
193  }
194  }
195  }
196  }
197  }
198 
199  // Do same for cyclics.
200  for (const polyPatch& pp : oldPatches)
201  {
202  if (isA<cyclicPolyPatch>(pp))
203  {
204  const cyclicPolyPatch& cycPatch =
205  refCast<const cyclicPolyPatch>(pp);
206 
207  if (!nCellsUsingFace.empty())
208  {
209  forAll(cycPatch, i)
210  {
211  label thisFacei = cycPatch.start() + i;
212  label otherFacei = cycPatch.transformGlobalFace(thisFacei);
213 
214  if
215  (
216  nCellsUsingFace[thisFacei] == 1
217  && nCellsUsingFace[otherFacei] == 0
218  )
219  {
220  nCellsUsingFace[thisFacei] = 3;
221  ++nUncoupled;
222  }
223  }
224  }
225  }
226  }
227 
228  if (syncPar)
229  {
230  reduce(nUncoupled, sumOp<label>());
231  }
232 
233  if (nUncoupled > 0)
234  {
235  Info<< "Uncoupled " << nUncoupled << " faces on coupled patches. "
236  << "(processorPolyPatch, cyclicPolyPatch)" << endl;
237  }
238 }
239 
240 
241 void Foam::fvMeshSubset::removeCellsImpl
242 (
243  const bitSet& cellsToRemove,
244  const labelList& exposedFaces,
245  const labelList& patchIDs,
246  const bool syncCouples
247 )
248 {
249  // Mesh changing engine.
250  polyTopoChange meshMod(baseMesh());
251 
252  removeCells cellRemover(baseMesh(), syncCouples);
253 
254  cellRemover.setRefinement
255  (
256  cellsToRemove,
257  exposedFaces,
258  patchIDs,
259  meshMod
260  );
261 
262  // Create mesh, return map from old to new mesh.
263  autoPtr<mapPolyMesh> map = meshMod.makeMesh
264  (
265  fvMeshSubsetPtr_,
266  IOobject
267  (
268  baseMesh().name(),
269  baseMesh().time().timeName(),
270  baseMesh().time(),
271  IOobject::READ_IF_PRESENT, // read fv* if present
273  ),
274  baseMesh(),
275  syncCouples
276  );
277 
278  pointMap_ = map().pointMap();
279  faceMap_ = map().faceMap();
280  cellMap_ = map().cellMap();
281  patchMap_ = identity(baseMesh().boundaryMesh().size());
282 }
283 
284 
285 Foam::labelList Foam::fvMeshSubset::subsetSubset
286 (
287  const label nElems,
288  const labelUList& selectedElements,
289  const labelUList& subsetMap
290 )
291 {
292  // Mark selected elements.
293  const bitSet selected(nElems, selectedElements);
294 
295  // Count subset of selected elements
296  label n = 0;
297  forAll(subsetMap, i)
298  {
299  if (selected[subsetMap[i]])
300  {
301  ++n;
302  }
303  }
304 
305  // Collect selected elements
306  labelList subsettedElements(n);
307  n = 0;
308 
309  forAll(subsetMap, i)
310  {
311  if (selected[subsetMap[i]])
312  {
313  subsettedElements[n] = i;
314  ++n;
315  }
316  }
317 
318  return subsettedElements;
319 }
320 
321 
322 void Foam::fvMeshSubset::subsetZones()
323 {
324  // Keep all zones, even if zero size.
325 
326  // PointZones
327 
328  const pointZoneMesh& pointZones = baseMesh().pointZones();
329 
330  List<pointZone*> pZonePtrs(pointZones.size());
331 
332  forAll(pointZones, zonei)
333  {
334  const pointZone& pz = pointZones[zonei];
335 
336  pZonePtrs[zonei] = new pointZone
337  (
338  pz.name(),
339  subsetSubset(baseMesh().nPoints(), pz, pointMap()),
340  zonei,
341  fvMeshSubsetPtr_().pointZones()
342  );
343  }
344 
345 
346  // FaceZones
347  // Do we need to remove zones where the side we're interested in
348  // no longer exists? Guess not.
349 
350  const faceZoneMesh& faceZones = baseMesh().faceZones();
351 
352  List<faceZone*> fZonePtrs(faceZones.size());
353 
354  forAll(faceZones, zonei)
355  {
356  const faceZone& fz = faceZones[zonei];
357 
358  // Expand faceZone to full mesh
359  // +1 : part of faceZone, flipped
360  // -1 : ,, , unflipped
361  // 0 : not part of faceZone
362  labelList zone(baseMesh().nFaces(), Zero);
363  forAll(fz, j)
364  {
365  if (fz.flipMap()[j])
366  {
367  zone[fz[j]] = 1;
368  }
369  else
370  {
371  zone[fz[j]] = -1;
372  }
373  }
374 
375  // Select faces
376  label nSub = 0;
377  forAll(faceMap(), j)
378  {
379  if (zone[faceMap()[j]] != 0)
380  {
381  ++nSub;
382  }
383  }
384  labelList subAddressing(nSub);
385  boolList subFlipStatus(nSub);
386  nSub = 0;
387  forAll(faceMap(), subFacei)
388  {
389  const label meshFacei = faceMap()[subFacei];
390  if (zone[meshFacei] != 0)
391  {
392  subAddressing[nSub] = subFacei;
393  const label subOwner = subMesh().faceOwner()[subFacei];
394  const label baseOwner = baseMesh().faceOwner()[meshFacei];
395  // If subowner is the same cell as the base keep the flip status
396  const bool sameOwner = (cellMap()[subOwner] == baseOwner);
397  const bool flip = (zone[meshFacei] == 1);
398  subFlipStatus[nSub] = (sameOwner == flip);
399 
400  ++nSub;
401  }
402  }
403 
404  fZonePtrs[zonei] = new faceZone
405  (
406  fz.name(),
407  subAddressing,
408  subFlipStatus,
409  zonei,
410  fvMeshSubsetPtr_().faceZones()
411  );
412  }
413 
414  // Cell Zones
415 
416  const cellZoneMesh& cellZones = baseMesh().cellZones();
417 
418  List<cellZone*> cZonePtrs(cellZones.size());
419 
420  forAll(cellZones, zonei)
421  {
422  const cellZone& cz = cellZones[zonei];
423 
424  cZonePtrs[zonei] = new cellZone
425  (
426  cz.name(),
427  subsetSubset(baseMesh().nCells(), cz, cellMap()),
428  zonei,
429  fvMeshSubsetPtr_().cellZones()
430  );
431  }
432 
433 
434  // Add the zones
435  fvMeshSubsetPtr_().addZones(pZonePtrs, fZonePtrs, cZonePtrs);
436 }
437 
438 
439 Foam::bitSet Foam::fvMeshSubset::getCellsToRemove
440 (
441  const bitSet& selectedCells
442 ) const
443 {
444  // Work on a copy
445  bitSet cellsToRemove(selectedCells);
446 
447  // Ensure we have the full range
448  cellsToRemove.resize(baseMesh().nCells(), false);
449 
450  // Invert the selection
451  cellsToRemove.flip();
452 
453  return cellsToRemove;
454 }
455 
456 
457 Foam::bitSet Foam::fvMeshSubset::getCellsToRemove
458 (
459  const label regioni,
460  const labelUList& regions
461 ) const
462 {
463  return BitSetOps::create
464  (
465  baseMesh().nCells(),
466  regioni,
467  regions,
468  false // on=false: invert return cells to remove
469  );
470 }
471 
472 
473 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
474 
475 Foam::fvMeshSubset::fvMeshSubset(const fvMesh& baseMesh)
476 :
477  baseMesh_(baseMesh),
478  fvMeshSubsetPtr_(nullptr),
479  faceFlipMapPtr_(nullptr),
480  pointMap_(),
481  faceMap_(),
482  cellMap_(),
483  patchMap_()
484 {}
485 
486 
487 Foam::fvMeshSubset::fvMeshSubset
488 (
489  const fvMesh& baseMesh,
490  const bitSet& selectedCells,
491  const label patchID,
492  const bool syncCouples
493 )
494 :
495  fvMeshSubset(baseMesh)
496 {
497  setCellSubset(selectedCells, patchID, syncCouples);
498 }
499 
500 
501 Foam::fvMeshSubset::fvMeshSubset
502 (
503  const fvMesh& baseMesh,
504  const labelHashSet& selectedCells,
505  const label patchID,
506  const bool syncCouples
507 )
508 :
509  fvMeshSubset(baseMesh)
510 {
511  setCellSubset(selectedCells, patchID, syncCouples);
512 }
513 
514 
515 Foam::fvMeshSubset::fvMeshSubset
516 (
517  const fvMesh& baseMesh,
518  const labelUList& selectedCells,
519  const label patchID,
520  const bool syncCouples
521 )
522 :
523  fvMeshSubset(baseMesh)
524 {
525  setCellSubset(selectedCells, patchID, syncCouples);
526 }
527 
528 
529 Foam::fvMeshSubset::fvMeshSubset
530 (
531  const fvMesh& baseMesh,
532  const label regioni,
533  const labelUList& regions,
534  const label patchID,
535  const bool syncCouples
536 )
537 :
538  fvMeshSubset(baseMesh)
539 {
540  setCellSubset(regioni, regions, patchID, syncCouples);
541 }
542 
543 
544 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
545 
547 {
548  fvMeshSubsetPtr_.reset(nullptr);
549  faceFlipMapPtr_.reset(nullptr);
550 
551  pointMap_.clear();
552  faceMap_.clear();
553  cellMap_.clear();
554  patchMap_.clear();
555 }
556 
557 
559 (
560  const bitSet& selectedCells,
561  const label patchID,
562  const bool syncPar
563 )
564 {
565  const cellList& oldCells = baseMesh().cells();
566  const faceList& oldFaces = baseMesh().faces();
567  const pointField& oldPoints = baseMesh().points();
568  const labelList& oldOwner = baseMesh().faceOwner();
569  const labelList& oldNeighbour = baseMesh().faceNeighbour();
570  const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
571  const label oldNInternalFaces = baseMesh().nInternalFaces();
572 
573  // Initial check on patches before doing anything time consuming.
574 
575  label wantedPatchID = patchID;
576 
577  if (wantedPatchID == -1)
578  {
579  // No explicit patch specified. Put in oldInternalFaces patch.
580  // Check if patch with this name already exists.
581  wantedPatchID = oldPatches.findPatchID(exposedPatchName);
582  }
583  else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
584  {
586  << "Non-existing patch index " << wantedPatchID << endl
587  << "Should be between 0 and " << oldPatches.size()-1
588  << abort(FatalError);
589  }
590 
591  // Clear all old maps and pointers
592  clear();
593 
594  // The selected cells - sorted in ascending order
595  cellMap_ = selectedCells.sortedToc();
596 
597  // The selectedCells should normally be in the [0,nCells) range,
598  // but it is better not to trust that.
599  {
600  label len = cellMap_.size();
601  for
602  (
603  label i = len-1;
604  i >= 0 && (cellMap_[i] >= oldCells.size());
605  --i
606  )
607  {
608  len = i;
609  }
610  cellMap_.resize(len);
611  }
612 
613 
614  // Mark all used faces. Count number of cells using them
615  // 0: face not used anymore
616  // 1: face used by one cell, face becomes/stays boundary face
617  // 2: face still used and remains internal face
618  // 3: face coupled and used by one cell only (so should become normal,
619  // non-coupled patch face)
620  //
621  // Note that this is not really necessary - but means we can size things
622  // correctly. Also makes handling coupled faces much easier.
623 
624  labelList nCellsUsingFace(oldFaces.size(), Zero);
625 
626  label nFacesInSet = 0;
627  forAll(oldFaces, oldFacei)
628  {
629  bool faceUsed = false;
630 
631  if (selectedCells.test(oldOwner[oldFacei]))
632  {
633  ++nCellsUsingFace[oldFacei];
634  faceUsed = true;
635  }
636 
637  if
638  (
639  baseMesh().isInternalFace(oldFacei)
640  && selectedCells.test(oldNeighbour[oldFacei])
641  )
642  {
643  ++nCellsUsingFace[oldFacei];
644  faceUsed = true;
645  }
646 
647  if (faceUsed)
648  {
649  ++nFacesInSet;
650  }
651  }
652  faceMap_.setSize(nFacesInSet);
653 
654  // Handle coupled faces. Modifies patch faces to be uncoupled to 3.
655  doCoupledPatches(syncPar, nCellsUsingFace);
656 
657 
658  // See which patch to use for exposed internal faces.
659  label oldInternalPatchID = 0;
660 
661  // Insert faces before which patch
662  label nextPatchID = oldPatches.size();
663 
664  // old to new patches
665  labelList globalPatchMap(oldPatches.size());
666 
667  // New patch size
668  label nbSize = oldPatches.size();
669 
670  if (wantedPatchID == -1)
671  {
672  // Create 'oldInternalFaces' patch at the end (or before
673  // processorPatches)
674  // and put all exposed internal faces in there.
675 
676  forAll(oldPatches, patchi)
677  {
678  if (isA<processorPolyPatch>(oldPatches[patchi]))
679  {
680  nextPatchID = patchi;
681  break;
682  }
683  ++oldInternalPatchID;
684  }
685 
686  ++nbSize;
687 
688  // adapt old to new patches for inserted patch
689  for (label oldPatchi = 0; oldPatchi < nextPatchID; oldPatchi++)
690  {
691  globalPatchMap[oldPatchi] = oldPatchi;
692  }
693  for
694  (
695  label oldPatchi = nextPatchID;
696  oldPatchi < oldPatches.size();
697  oldPatchi++
698  )
699  {
700  globalPatchMap[oldPatchi] = oldPatchi+1;
701  }
702  }
703  else
704  {
705  oldInternalPatchID = wantedPatchID;
706  nextPatchID = wantedPatchID+1;
707 
708  // old to new patches
709  globalPatchMap = identity(oldPatches.size());
710  }
711 
712  labelList boundaryPatchSizes(nbSize, Zero);
713 
714 
715  // Make a global-to-local point map
716  labelList globalPointMap(oldPoints.size(), -1);
717  labelList globalFaceMap(oldFaces.size(), -1);
718 
719  label facei = 0;
720 
721  // 1. Pick up all preserved internal faces.
722  for (label oldFacei = 0; oldFacei < oldNInternalFaces; ++oldFacei)
723  {
724  if (nCellsUsingFace[oldFacei] == 2)
725  {
726  globalFaceMap[oldFacei] = facei;
727  faceMap_[facei++] = oldFacei;
728 
729  // Mark all points from the face
730  markUsed(oldFaces[oldFacei], globalPointMap);
731  }
732  }
733 
734  // These are all the internal faces in the mesh.
735  const label nInternalFaces = facei;
736 
737  // 2. Boundary faces up to where we want to insert old internal faces
738  for
739  (
740  label oldPatchi = 0;
741  oldPatchi < oldPatches.size()
742  && oldPatchi < nextPatchID;
743  oldPatchi++
744  )
745  {
746  const polyPatch& oldPatch = oldPatches[oldPatchi];
747 
748  label oldFacei = oldPatch.start();
749 
750  forAll(oldPatch, i)
751  {
752  if (nCellsUsingFace[oldFacei] == 1)
753  {
754  // Boundary face is kept.
755 
756  // Mark face and increment number of points in set
757  globalFaceMap[oldFacei] = facei;
758  faceMap_[facei++] = oldFacei;
759 
760  // Mark all points from the face
761  markUsed(oldFaces[oldFacei], globalPointMap);
762 
763  // Increment number of patch faces
764  ++boundaryPatchSizes[globalPatchMap[oldPatchi]];
765  }
766  ++oldFacei;
767  }
768  }
769 
770  // 3a. old internal faces that have become exposed.
771  for (label oldFacei = 0; oldFacei < oldNInternalFaces; ++oldFacei)
772  {
773  if (nCellsUsingFace[oldFacei] == 1)
774  {
775  globalFaceMap[oldFacei] = facei;
776  faceMap_[facei++] = oldFacei;
777 
778  // Mark all points from the face
779  markUsed(oldFaces[oldFacei], globalPointMap);
780 
781  // Increment number of patch faces
782  ++boundaryPatchSizes[oldInternalPatchID];
783  }
784  }
785 
786  // 3b. coupled patch faces that have become uncoupled.
787  for
788  (
789  label oldFacei = oldNInternalFaces;
790  oldFacei < oldFaces.size();
791  oldFacei++
792  )
793  {
794  if (nCellsUsingFace[oldFacei] == 3)
795  {
796  globalFaceMap[oldFacei] = facei;
797  faceMap_[facei++] = oldFacei;
798 
799  // Mark all points from the face
800  markUsed(oldFaces[oldFacei], globalPointMap);
801 
802  // Increment number of patch faces
803  ++boundaryPatchSizes[oldInternalPatchID];
804  }
805  }
806 
807  // 4. Remaining boundary faces
808  for
809  (
810  label oldPatchi = nextPatchID;
811  oldPatchi < oldPatches.size();
812  oldPatchi++
813  )
814  {
815  const polyPatch& oldPatch = oldPatches[oldPatchi];
816 
817  label oldFacei = oldPatch.start();
818 
819  forAll(oldPatch, i)
820  {
821  if (nCellsUsingFace[oldFacei] == 1)
822  {
823  // Boundary face is kept.
824 
825  // Mark face and increment number of points in set
826  globalFaceMap[oldFacei] = facei;
827  faceMap_[facei++] = oldFacei;
828 
829  // Mark all points from the face
830  markUsed(oldFaces[oldFacei], globalPointMap);
831 
832  // Increment number of patch faces
833  ++boundaryPatchSizes[globalPatchMap[oldPatchi]];
834  }
835  ++oldFacei;
836  }
837  }
838 
839  if (facei != nFacesInSet)
840  {
842  << "Problem" << abort(FatalError);
843  }
844 
845 
846  // Grab the points map
847  label nPointsInSet = 0;
848 
849  forAll(globalPointMap, pointi)
850  {
851  if (globalPointMap[pointi] != -1)
852  {
853  ++nPointsInSet;
854  }
855  }
856  pointMap_.setSize(nPointsInSet);
857 
858  nPointsInSet = 0;
859 
860  forAll(globalPointMap, pointi)
861  {
862  if (globalPointMap[pointi] != -1)
863  {
864  pointMap_[nPointsInSet] = pointi;
865  globalPointMap[pointi] = nPointsInSet;
866  ++nPointsInSet;
867  }
868  }
869 
870 
871  //Pout<< "Number of points,faces,cells in new mesh : "
872  // << pointMap_.size() << ' '
873  // << faceMap_.size() << ' '
874  // << cellMap_.size() << nl;
875 
876 
877  // Make a new mesh
878 
879  //
880  // Create new points
881  //
882  pointField newPoints(pointUIndList(oldPoints, pointMap_));
883 
884 
885  //
886  // Create new faces
887  //
888 
889  faceList newFaces(faceMap_.size());
890  {
891  auto iter = newFaces.begin();
892  const auto& renumbering = globalPointMap;
893 
894  // Internal faces
895  for (label facei = 0; facei < nInternalFaces; ++facei)
896  {
897  face& newItem = *iter;
898  ++iter;
899 
900  const face& oldItem = oldFaces[faceMap_[facei]];
901 
902  // Copy relabelled
903  newItem.resize(oldItem.size());
904 
905  forAll(oldItem, i)
906  {
907  newItem[i] = renumbering[oldItem[i]];
908  }
909  }
910 
911  // Boundary faces - may need to be flipped
912  for (label facei = nInternalFaces; facei < faceMap_.size(); ++facei)
913  {
914  const label oldFacei = faceMap_[facei];
915 
916  face& newItem = *iter;
917  ++iter;
918 
919  const face oldItem =
920  (
921  (
922  baseMesh().isInternalFace(oldFacei)
923  && selectedCells.test(oldNeighbour[oldFacei])
924  && !selectedCells.test(oldOwner[oldFacei])
925  )
926  // Face flipped to point outwards
927  ? oldFaces[oldFacei].reverseFace()
928  : oldFaces[oldFacei]
929  );
930 
931  // Copy relabelled
932  newItem.resize(oldItem.size());
933 
934  forAll(oldItem, i)
935  {
936  newItem[i] = renumbering[oldItem[i]];
937  }
938  }
939  }
940 
941 
942  //
943  // Create new cells
944  //
945 
946  cellList newCells(cellMap_.size());
947  {
948  auto iter = newCells.begin();
949  const auto& renumbering = globalFaceMap;
950 
951  for (const label oldCelli : cellMap_)
952  {
953  cell& newItem = *iter;
954  ++iter;
955 
956  const labelList& oldItem = oldCells[oldCelli];
957 
958  // Copy relabelled
959  newItem.resize(oldItem.size());
960 
961  forAll(oldItem, i)
962  {
963  newItem[i] = renumbering[oldItem[i]];
964  }
965  }
966  }
967 
968 
969  // Make a new mesh
970  //
971  // Note that mesh gets registered with same name as original mesh.
972  // This is not proper but cannot be avoided since otherwise
973  // surfaceInterpolation cannot find its fvSchemes.
974  // It will try to read, for example. "system/region0SubSet/fvSchemes"
975  //
976  fvMeshSubsetPtr_ = autoPtr<fvMesh>::New
977  (
978  IOobject
979  (
980  baseMesh().name(),
981  baseMesh().time().timeName(),
982  baseMesh().time(),
983  IOobject::NO_READ, // do not read any dictionaries
985  ),
986  baseMesh(), // get dictionaries from base mesh
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  if (oldInternalPatchID != oldPatchi)
1060  {
1061  // Pure subset of patch faces (no internal faces added to this
1062  // patch). Can use mapping.
1063  labelList map(newSize);
1064  for (label patchFacei = 0; patchFacei < newSize; patchFacei++)
1065  {
1066  const label facei = patchStart+patchFacei;
1067  const label oldFacei = faceMap_[facei];
1068  map[patchFacei] = oldPatches[oldPatchi].whichFace(oldFacei);
1069  }
1070 
1071  newBoundary[nNewPatches] = oldPatches[oldPatchi].clone
1072  (
1073  fvMeshSubsetPtr_().boundaryMesh(),
1074  nNewPatches,
1075  map,
1076  patchStart
1077  ).ptr();
1078  }
1079  else
1080  {
1081  // Clone (even if 0 size)
1082  newBoundary[nNewPatches] = oldPatches[oldPatchi].clone
1083  (
1084  fvMeshSubsetPtr_().boundaryMesh(),
1085  nNewPatches,
1086  newSize,
1087  patchStart
1088  ).ptr();
1089  }
1090 
1091  patchStart += newSize;
1092  patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1093  ++nNewPatches;
1094  }
1095 
1096  // Inserted patch
1097 
1098  if (wantedPatchID == -1)
1099  {
1100  label oldInternalSize = boundaryPatchSizes[oldInternalPatchID];
1101 
1102  if (syncPar)
1103  {
1104  reduce(oldInternalSize, sumOp<label>());
1105  }
1106 
1107  // Newly created patch so is at end. Check if any faces in it.
1108  if (oldInternalSize > 0)
1109  {
1110  newBoundary[nNewPatches] = new emptyPolyPatch
1111  (
1112  exposedPatchName,
1113  boundaryPatchSizes[oldInternalPatchID],
1114  patchStart,
1115  nNewPatches,
1116  fvMeshSubsetPtr_().boundaryMesh(),
1117  emptyPolyPatch::typeName
1118  );
1119 
1120  //Pout<< " " << exposedPatchName << " : "
1121  // << boundaryPatchSizes[oldInternalPatchID] << endl;
1122 
1123  // The index for the first patch is -1 as it originates from
1124  // the internal faces
1125  patchStart += boundaryPatchSizes[oldInternalPatchID];
1126  patchMap_[nNewPatches] = -1;
1127  ++nNewPatches;
1128  }
1129  }
1130 
1131  // Old patches
1132 
1133  for
1134  (
1135  label oldPatchi = nextPatchID;
1136  oldPatchi < oldPatches.size();
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[nNewPatches] = oldPatches[oldPatchi].clone
1155  (
1156  fvMeshSubsetPtr_().boundaryMesh(),
1157  nNewPatches,
1158  map,
1159  patchStart
1160  ).ptr();
1161  }
1162  else
1163  {
1164  // Patch still exists. Add it
1165  newBoundary[nNewPatches] = oldPatches[oldPatchi].clone
1166  (
1167  fvMeshSubsetPtr_().boundaryMesh(),
1168  nNewPatches,
1169  newSize,
1170  patchStart
1171  ).ptr();
1172  }
1173 
1174  //Pout<< " " << oldPatches[oldPatchi].name() << " : "
1175  // << newSize << endl;
1176 
1177  patchStart += newSize;
1178  patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1179  ++nNewPatches;
1180  }
1181 
1182 
1183  // Reset the patch lists
1184  newBoundary.setSize(nNewPatches);
1185  patchMap_.setSize(nNewPatches);
1186 
1187 
1188  // Add the fvPatches
1189  fvMeshSubsetPtr_().addFvPatches(newBoundary, syncPar);
1190 
1191  // Subset and add any zones
1192  subsetZones();
1193 }
1194 
1195 
1198  const labelUList& selectedCells,
1199  const label patchID,
1200  const bool syncPar
1201 )
1202 {
1203  setCellSubset
1204  (
1205  BitSetOps::create(baseMesh().nCells(), selectedCells),
1206  patchID,
1207  syncPar
1208  );
1209 }
1210 
1211 
1214  const labelHashSet& selectedCells,
1215  const label patchID,
1216  const bool syncPar
1217 )
1218 {
1219  setCellSubset
1220  (
1221  BitSetOps::create(baseMesh().nCells(), selectedCells),
1222  patchID,
1223  syncPar
1224  );
1225 }
1226 
1227 
1230  const label regioni,
1231  const labelUList& regions,
1232  const label patchID,
1233  const bool syncPar
1234 )
1235 {
1236  setCellSubset
1237  (
1238  BitSetOps::create(baseMesh().nCells(), regioni, regions),
1239  patchID,
1240  syncPar
1241  );
1242 }
1243 
1244 
1247  const bitSet& selectedCells,
1248  const bool syncCouples
1249 ) const
1250 {
1251  return
1252  removeCells(baseMesh(), syncCouples)
1253  .getExposedFaces(getCellsToRemove(selectedCells));
1254 }
1255 
1256 
1259  const label regioni,
1260  const labelUList& regions,
1261  const bool syncCouples
1262 ) const
1263 {
1264  return
1265  removeCells(baseMesh(), syncCouples)
1266  .getExposedFaces(getCellsToRemove(regioni, regions));
1267 }
1268 
1269 
1272  const bitSet& selectedCells,
1273  const labelList& exposedFaces,
1274  const labelList& patchIDs,
1275  const bool syncCouples
1276 )
1277 {
1278  removeCellsImpl
1279  (
1280  getCellsToRemove(selectedCells),
1281  exposedFaces,
1282  patchIDs,
1283  syncCouples
1284  );
1285 }
1286 
1287 
1290  const label selectRegion,
1291  const labelList& regions,
1292  const labelList& exposedFaces,
1293  const labelList& patchIDs,
1294  const bool syncCouples
1295 )
1296 {
1297  removeCellsImpl
1298  (
1299  getCellsToRemove(selectRegion, regions),
1300  exposedFaces,
1301  patchIDs,
1302  syncCouples
1303  );
1304 }
1305 
1306 
1307 // ************************************************************************* //
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:71
BitOps.H
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::BitSetOps::create
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:123
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: blockMeshMergeTopological.C:94
boolList.H
Foam::labelMax
constexpr label labelMax
Definition: label.H:61
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:63
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 (0)
Definition: zero.H:131
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
polyTopoChange.H
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
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:492
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:559
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:69
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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:520
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:546
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:486
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::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:68
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:480
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:1107
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
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:83
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:313
emptyPolyPatch.H
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:1094
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
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:1246
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:464
clear
patchWriters clear()
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
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:62
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
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:80
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:532
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:446
Foam::removeCells::getExposedFaces
labelList getExposedFaces(const bitSet &removedCell) const
Get labels of faces exposed after cells removal.
Definition: removeCells.C:97