fvMeshSubset.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  Post-processing mesh subset tool. Given the original mesh and the
26  list of selected cells, it creates the mesh consisting only of the
27  desired cells, with the mapping list for points, faces, and cells.
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "fvMeshSubset.H"
32 #include "boolList.H"
33 #include "Pstream.H"
34 #include "internalPolyPatch.H"
35 #include "demandDrivenData.H"
36 #include "cyclicPolyPatch.H"
37 #include "processorPolyPatch.H"
38 #include "removeCells.H"
39 #include "polyTopoChange.H"
40 #include "polyTopoChangeMap.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 bool Foam::fvMeshSubset::checkCellSubset() const
53 {
54  if (fvMeshSubsetPtr_.empty())
55  {
57  << "void setCellSubset(const labelHashSet& cellsToSubset)" << endl
58  << "before attempting to access subset data"
59  << abort(FatalError);
60 
61  return false;
62  }
63  else
64  {
65  return true;
66  }
67 }
68 
69 
70 void Foam::fvMeshSubset::markPoints
71 (
72  const labelList& curPoints,
73  Map<label>& pointMap
74 )
75 {
76  forAll(curPoints, pointi)
77  {
78  // Note: insert will only insert if not yet there.
79  pointMap.insert(curPoints[pointi], 0);
80  }
81 }
82 
83 
84 void Foam::fvMeshSubset::markPoints
85 (
86  const labelList& curPoints,
87  labelList& pointMap
88 )
89 {
90  forAll(curPoints, pointi)
91  {
92  pointMap[curPoints[pointi]] = 0;
93  }
94 }
95 
96 
97 void Foam::fvMeshSubset::doCoupledPatches
98 (
99  const bool syncPar,
100  Map<label>& facesToSubset,
101  labelList& nCellsUsingFace
102 ) const
103 {
104  // Synchronise facesToSubset on both sides of coupled patches.
105  // Marks faces that become 'uncoupled' with 3.
106 
107  const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
108 
109  label nUncoupled = 0;
110 
111  if (syncPar && Pstream::parRun())
112  {
113  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
114 
115  // Send face usage across processor patches
116  forAll(oldPatches, oldPatchi)
117  {
118  const polyPatch& pp = oldPatches[oldPatchi];
119 
120  if (isA<processorPolyPatch>(pp))
121  {
122  const processorPolyPatch& procPatch =
123  refCast<const processorPolyPatch>(pp);
124 
125  UOPstream toNeighbour(procPatch.neighbProcNo(), pBufs);
126 
127  if (!facesToSubset.empty())
128  {
129  DynamicList<label> patchFacesToSubset;
130  forAll(pp, i)
131  {
132  if
133  (
134  facesToSubset.found(pp.start()+i)
135  && facesToSubset[pp.start()+i] == 1
136  )
137  {
138  patchFacesToSubset.append(i);
139  }
140  }
141  toNeighbour << patchFacesToSubset;
142  }
143  else if (!nCellsUsingFace.empty())
144  {
145  toNeighbour <<
146  SubList<label>(nCellsUsingFace, pp.size(), pp.start());
147  }
148  else
149  {
150  toNeighbour << labelList();
151  }
152  }
153  }
154 
155  pBufs.finishedSends();
156 
157  // Receive face usage count and check for faces that become uncoupled.
158  forAll(oldPatches, oldPatchi)
159  {
160  const polyPatch& pp = oldPatches[oldPatchi];
161 
162  if (isA<processorPolyPatch>(pp))
163  {
164  const processorPolyPatch& procPatch =
165  refCast<const processorPolyPatch>(pp);
166 
167  UIPstream fromNeighbour(procPatch.neighbProcNo(), pBufs);
168 
169  const labelList nbrList(fromNeighbour);
170 
171  // Combine with this side.
172 
173  if (!facesToSubset.empty())
174  {
175  const labelHashSet nbrPatchFacesToSubset(nbrList);
176 
177  forAll(pp, i)
178  {
179  if
180  (
181  facesToSubset.found(pp.start()+i)
182  && facesToSubset[pp.start()+i] == 1
183  && !nbrPatchFacesToSubset.found(i)
184  )
185  {
186  // Face's neighbour is no longer there. Mark face
187  // off as coupled
188  facesToSubset[pp.start()+i] = 3;
189  nUncoupled++;
190  }
191  }
192  }
193  else if (!nCellsUsingFace.empty())
194  {
195  const labelList& nbrCellsUsingFace(nbrList);
196 
197  // Combine with this side.
198 
199  forAll(pp, i)
200  {
201  if
202  (
203  nCellsUsingFace[pp.start()+i] == 1
204  && nbrCellsUsingFace[i] == 0
205  )
206  {
207  // Face's neighbour is no longer there. Mark face
208  // off as coupled
209  nCellsUsingFace[pp.start()+i] = 3;
210  nUncoupled++;
211  }
212  }
213  }
214  }
215  }
216  }
217 
218  // Do same for cyclics.
219  forAll(oldPatches, oldPatchi)
220  {
221  const polyPatch& pp = oldPatches[oldPatchi];
222 
223  if (isA<cyclicPolyPatch>(pp))
224  {
225  const cyclicPolyPatch& cycPatch =
226  refCast<const cyclicPolyPatch>(pp);
227 
228  if (!facesToSubset.empty())
229  {
230  forAll(cycPatch, i)
231  {
232  const label thisFacei = cycPatch.start() + i;
233  const label otherFacei =
234  cycPatch.transformGlobalFace(thisFacei);
235 
236  if
237  (
238  facesToSubset.found(thisFacei)
239  && facesToSubset[thisFacei] == 1
240  && !facesToSubset.found(otherFacei)
241  )
242  {
243  facesToSubset[thisFacei] = 3;
244  nUncoupled++;
245  }
246  }
247  }
248  else if (!nCellsUsingFace.empty())
249  {
250  forAll(cycPatch, i)
251  {
252  const label thisFacei = cycPatch.start() + i;
253  const label otherFacei =
254  cycPatch.transformGlobalFace(thisFacei);
255 
256  if
257  (
258  nCellsUsingFace[thisFacei] == 1
259  && nCellsUsingFace[otherFacei] == 0
260  )
261  {
262  nCellsUsingFace[thisFacei] = 3;
263  nUncoupled++;
264  }
265  }
266  }
267  }
268  }
269 
270  if (syncPar)
271  {
272  reduce(nUncoupled, sumOp<label>());
273  }
274 
275  if (debug && nUncoupled > 0)
276  {
277  Info<< "Uncoupled " << nUncoupled << " faces on coupled patches. "
278  << "(processorPolyPatch, cyclicPolyPatch)" << endl;
279  }
280 }
281 
282 
283 Foam::labelList Foam::fvMeshSubset::subset
284 (
285  const label nElems,
286  const labelList& selectedElements,
287  const labelList& subsetMap
288 )
289 {
290  // Mark selected elements.
291  boolList selected(nElems, false);
292  forAll(selectedElements, i)
293  {
294  selected[selectedElements[i]] = true;
295  }
296 
297  // Count subset of selected elements
298  label n = 0;
299  forAll(subsetMap, i)
300  {
301  if (selected[subsetMap[i]])
302  {
303  n++;
304  }
305  }
306 
307  // Collect selected elements
308  labelList subsettedElements(n);
309  n = 0;
310 
311  forAll(subsetMap, i)
312  {
313  if (selected[subsetMap[i]])
314  {
315  subsettedElements[n++] = i;
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  const pointZoneList& pointZones = baseMesh().pointZones();
328 
329  // PointZones
330  List<pointZone*> pZonePtrs(pointZones.size());
331 
332  forAll(pointZones, i)
333  {
334  const pointZone& pz = pointZones[i];
335 
336  pZonePtrs[i] = new pointZone
337  (
338  pz.name(),
339  subset(baseMesh().nPoints(), pz, pointMap()),
340  fvMeshSubsetPtr_().pointZones()
341  );
342  }
343 
344 
345  // FaceZones
346  const faceZoneList& faceZones = baseMesh().faceZones();
347 
348  // Do we need to remove zones where the side we're interested in
349  // no longer exists? Guess not.
350  List<faceZone*> fZonePtrs(faceZones.size());
351 
352  forAll(faceZones, i)
353  {
354  const faceZone& fz = faceZones[i];
355 
356  // Expand faceZone to full mesh
357  // +1 : part of faceZone, flipped
358  // -1 : ,, , unflipped
359  // 0 : not part of faceZone
360  labelList zone(baseMesh().nFaces(), 0);
361  forAll(fz, j)
362  {
363  if (fz.flipMap()[j])
364  {
365  zone[fz[j]] = 1;
366  }
367  else
368  {
369  zone[fz[j]] = -1;
370  }
371  }
372 
373  // Select faces
374  label nSub = 0;
375  forAll(faceMap(), j)
376  {
377  if (zone[faceMap()[j]] != 0)
378  {
379  nSub++;
380  }
381  }
382  labelList subAddressing(nSub);
383  boolList subFlipStatus(nSub);
384  nSub = 0;
385  forAll(faceMap(), subFacei)
386  {
387  const label meshFacei = faceMap()[subFacei];
388  if (zone[meshFacei] != 0)
389  {
390  subAddressing[nSub] = subFacei;
391  const label subOwner = subMesh().faceOwner()[subFacei];
392  const label baseOwner = baseMesh().faceOwner()[meshFacei];
393  // If subowner is the same cell as the base keep the flip status
394  const bool sameOwner = (cellMap()[subOwner] == baseOwner);
395  const bool flip = (zone[meshFacei] == 1);
396  subFlipStatus[nSub] = (sameOwner == flip);
397 
398  nSub++;
399  }
400  }
401 
402  fZonePtrs[i] = new faceZone
403  (
404  fz.name(),
405  subAddressing,
406  subFlipStatus,
407  fvMeshSubsetPtr_().faceZones()
408  );
409  }
410 
411  // CellZones
412  const cellZoneList& cellZones = baseMesh().cellZones();
413 
414  List<cellZone*> cZonePtrs(cellZones.size());
415 
416  forAll(cellZones, i)
417  {
418  const cellZone& cz = cellZones[i];
419 
420  cZonePtrs[i] = new cellZone
421  (
422  cz.name(),
423  subset(baseMesh().nCells(), cz, cellMap()),
424  fvMeshSubsetPtr_().cellZones()
425  );
426  }
427 
428 
429  // Add the zones
430  fvMeshSubsetPtr_().addZones(pZonePtrs, fZonePtrs, cZonePtrs);
431 }
432 
433 
434 Foam::labelList Foam::fvMeshSubset::getCellsToRemove
435 (
436  const labelList& region,
437  const label currentRegion
438 ) const
439 {
440  // Count
441  label nKeep = 0;
442  forAll(region, cellI)
443  {
444  if (region[cellI] == currentRegion)
445  {
446  nKeep++;
447  }
448  }
449 
450  // Collect cells to remove
451  label nRemove = baseMesh().nCells() - nKeep;
452  labelList cellsToRemove(nRemove);
453 
454  nRemove = 0;
455  forAll(region, cellI)
456  {
457  if (region[cellI] != currentRegion)
458  {
459  cellsToRemove[nRemove++] = cellI;
460  }
461  }
462 
463  return cellsToRemove;
464 }
465 
466 
467 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
468 
470 :
471  baseMesh_(baseMesh),
472  fvMeshSubsetPtr_(nullptr),
473  pointMap_(0),
474  faceMap_(0),
475  cellMap_(0),
476  patchMap_(0),
477  faceFlipMapPtr_()
478 {}
479 
480 
481 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
482 
484 {}
485 
486 
487 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
488 
490 (
491  const labelHashSet& globalCellMap,
492  const label patchID,
493  const bool syncPar
494 )
495 {
496  // Initial check on patches before doing anything time consuming.
497  const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
498  const cellList& oldCells = baseMesh().cells();
499  const faceList& oldFaces = baseMesh().faces();
500  const pointField& oldPoints = baseMesh().points();
501  const labelList& oldOwner = baseMesh().faceOwner();
502  const labelList& oldNeighbour = baseMesh().faceNeighbour();
503 
504  label wantedPatchID = patchID;
505 
506  if (wantedPatchID == -1)
507  {
508  // No explicit patch specified. Put in oldInternalFaces patch.
509  // Check if patch with this name already exists.
510  wantedPatchID = oldPatches.findIndex("oldInternalFaces");
511  }
512  else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
513  {
515  << "Non-existing patch index " << wantedPatchID << endl
516  << "Should be between 0 and " << oldPatches.size()-1
517  << abort(FatalError);
518  }
519 
520 
521  // Clear demand driven data
522  faceFlipMapPtr_.clear();
523 
524 
525  cellMap_ = globalCellMap.toc();
526 
527  // Sort the cell map in the ascending order
528  sort(cellMap_);
529 
530  // Approximate sizing parameters for face and point lists
531  const label avgNFacesPerCell = 6;
532  const label avgNPointsPerFace = 4;
533 
534 
535  const label nCellsInSet = cellMap_.size();
536 
537  // Mark all used faces
538 
539  Map<label> facesToSubset(avgNFacesPerCell*nCellsInSet);
540 
541  forAll(cellMap_, celli)
542  {
543  // Mark all faces from the cell
544  const labelList& curFaces = oldCells[cellMap_[celli]];
545 
546  forAll(curFaces, facei)
547  {
548  if (!facesToSubset.found(curFaces[facei]))
549  {
550  facesToSubset.insert(curFaces[facei], 1);
551  }
552  else
553  {
554  facesToSubset[curFaces[facei]]++;
555  }
556  }
557  }
558 
559  // Handle coupled faces. Modifies patch faces to be uncoupled to 3.
560  labelList empty;
561  doCoupledPatches(syncPar, facesToSubset, empty);
562 
563  // Mark all used points and make a global-to-local face map
564  Map<label> globalFaceMap(facesToSubset.size());
565 
566  // Make a global-to-local point map
567  Map<label> globalPointMap(avgNPointsPerFace*facesToSubset.size());
568 
569  // This is done in two goes, so that the boundary faces are last
570  // in the list. Because of this, I need to create the face map
571  // along the way rather than just grab the table of contents.
572  labelList facesToc = facesToSubset.toc();
573  sort(facesToc);
574  faceMap_.setSize(facesToc.size());
575 
576  // 1. Get all faces that will be internal to the submesh.
577  forAll(facesToc, facei)
578  {
579  if (facesToSubset[facesToc[facei]] == 2)
580  {
581  // Mark face and increment number of points in set
582  faceMap_[globalFaceMap.size()] = facesToc[facei];
583  globalFaceMap.insert(facesToc[facei], globalFaceMap.size());
584 
585  // Mark all points from the face
586  markPoints(oldFaces[facesToc[facei]], globalPointMap);
587  }
588  }
589 
590  // These are all the internal faces in the mesh.
591  const label nInternalFaces = globalFaceMap.size();
592 
593 
594  // Where to insert old internal faces.
595  label oldPatchStart = labelMax;
596  if (wantedPatchID != -1)
597  {
598  oldPatchStart = oldPatches[wantedPatchID].start();
599  }
600 
601 
602  label facei = 0;
603 
604  // 2. Boundary faces up to where we want to insert old internal faces
605  for (; facei< facesToc.size(); facei++)
606  {
607  if (facesToc[facei] >= oldPatchStart)
608  {
609  break;
610  }
611  if
612  (
613  !baseMesh().isInternalFace(facesToc[facei])
614  && facesToSubset[facesToc[facei]] == 1
615  )
616  {
617  // Mark face and increment number of points in set
618  faceMap_[globalFaceMap.size()] = facesToc[facei];
619  globalFaceMap.insert(facesToc[facei], globalFaceMap.size());
620 
621  // Mark all points from the face
622  markPoints(oldFaces[facesToc[facei]], globalPointMap);
623  }
624  }
625 
626  // 3. old internal faces and uncoupled faces
627  forAll(facesToc, intFacei)
628  {
629  if
630  (
631  (
632  baseMesh().isInternalFace(facesToc[intFacei])
633  && facesToSubset[facesToc[intFacei]] == 1
634  )
635  || (
636  !baseMesh().isInternalFace(facesToc[intFacei])
637  && facesToSubset[facesToc[intFacei]] == 3
638  )
639  )
640  {
641  // Mark face and increment number of points in set
642  faceMap_[globalFaceMap.size()] = facesToc[intFacei];
643  globalFaceMap.insert(facesToc[intFacei], globalFaceMap.size());
644 
645  // Mark all points from the face
646  markPoints(oldFaces[facesToc[intFacei]], globalPointMap);
647  }
648  }
649 
650  // 4. Remaining boundary faces
651  for (; facei< facesToc.size(); facei++)
652  {
653  if
654  (
655  !baseMesh().isInternalFace(facesToc[facei])
656  && facesToSubset[facesToc[facei]] == 1
657  )
658  {
659  // Mark face and increment number of points in set
660  faceMap_[globalFaceMap.size()] = facesToc[facei];
661  globalFaceMap.insert(facesToc[facei], globalFaceMap.size());
662 
663  // Mark all points from the face
664  markPoints(oldFaces[facesToc[facei]], globalPointMap);
665  }
666  }
667 
668 
669 
670  // Grab the points map
671  pointMap_ = globalPointMap.toc();
672  sort(pointMap_);
673 
674  forAll(pointMap_, pointi)
675  {
676  globalPointMap[pointMap_[pointi]] = pointi;
677  }
678 
679  // Make a new mesh
680  pointField newPoints(globalPointMap.size());
681 
682  label nNewPoints = 0;
683 
684  forAll(pointMap_, pointi)
685  {
686  newPoints[nNewPoints] = oldPoints[pointMap_[pointi]];
687  nNewPoints++;
688  }
689 
690  faceList newFaces(globalFaceMap.size());
691 
692  label nNewFaces = 0;
693 
694  // Make internal faces
695  for (label facei = 0; facei < nInternalFaces; facei++)
696  {
697  const face& oldF = oldFaces[faceMap_[facei]];
698 
699  face newF(oldF.size());
700 
701  forAll(newF, i)
702  {
703  newF[i] = globalPointMap[oldF[i]];
704  }
705 
706  newFaces[nNewFaces] = newF;
707  nNewFaces++;
708  }
709 
710  // Make boundary faces
711 
712  label nbSize = oldPatches.size();
713  label oldInternalPatchID = -1;
714 
715  if (wantedPatchID == -1)
716  {
717  // Create 'oldInternalFaces' patch at the end
718  // and put all exposed internal faces in there.
719  oldInternalPatchID = nbSize;
720  nbSize++;
721 
722  }
723  else
724  {
725  oldInternalPatchID = wantedPatchID;
726  }
727 
728 
729  // Grad size and start of each patch on the fly. Because of the
730  // structure of the underlying mesh, the patches will appear in the
731  // ascending order
732  labelList boundaryPatchSizes(nbSize, 0);
733 
734  // Assign boundary faces. Visited in order of faceMap_.
735  for (label facei = nInternalFaces; facei < faceMap_.size(); facei++)
736  {
737  const label oldFacei = faceMap_[facei];
738 
739  face oldF = oldFaces[oldFacei];
740 
741  // Turn the faces as necessary to point outwards
742  if (baseMesh().isInternalFace(oldFacei))
743  {
744  // Internal face. Possibly turned the wrong way round
745  if
746  (
747  !globalCellMap.found(oldOwner[oldFacei])
748  && globalCellMap.found(oldNeighbour[oldFacei])
749  )
750  {
751  oldF = oldFaces[oldFacei].reverseFace();
752  }
753 
754  // Update count for patch
755  boundaryPatchSizes[oldInternalPatchID]++;
756  }
757  else if (facesToSubset[oldFacei] == 3)
758  {
759  // Uncoupled face. Increment the old patch.
760  boundaryPatchSizes[oldInternalPatchID]++;
761  }
762  else
763  {
764  // Boundary face. Increment the appropriate patch
765  const label patchOfFace = oldPatches.whichPatch(oldFacei);
766 
767  // Update count for patch
768  boundaryPatchSizes[patchOfFace]++;
769  }
770 
771  face newF(oldF.size());
772 
773  forAll(newF, i)
774  {
775  newF[i] = globalPointMap[oldF[i]];
776  }
777 
778  newFaces[nNewFaces] = newF;
779  nNewFaces++;
780  }
781 
782 
783 
784  // Create cells
785  cellList newCells(nCellsInSet);
786 
787  label nNewCells = 0;
788 
789  forAll(cellMap_, celli)
790  {
791  const labelList& oldC = oldCells[cellMap_[celli]];
792 
793  labelList newC(oldC.size());
794 
795  forAll(newC, i)
796  {
797  newC[i] = globalFaceMap[oldC[i]];
798  }
799 
800  newCells[nNewCells] = cell(newC);
801  nNewCells++;
802  }
803 
804 
805  // Delete any old mesh
806  fvMeshSubsetPtr_.clear();
807 
808  // Make a new mesh
809  fvMeshSubsetPtr_.reset
810  (
811  new fvMesh
812  (
813  IOobject
814  (
815  baseMesh().name() + "SubSet",
816  baseMesh().time().name(),
817  baseMesh().time(),
820  ),
821  move(newPoints),
822  move(newFaces),
823  move(newCells)
824  )
825  );
826 
827 
828  // Add old patches
829  List<polyPatch*> newBoundary(nbSize);
830  patchMap_.setSize(nbSize);
831  label nNewPatches = 0;
832  label patchStart = nInternalFaces;
833 
834 
835  forAll(oldPatches, patchi)
836  {
837  if (boundaryPatchSizes[patchi] > 0)
838  {
839  // Patch still exists. Add it
840  newBoundary[nNewPatches] = oldPatches[patchi].clone
841  (
842  fvMeshSubsetPtr_().boundaryMesh(),
843  nNewPatches,
844  boundaryPatchSizes[patchi],
845  patchStart
846  ).ptr();
847 
848  patchStart += boundaryPatchSizes[patchi];
849  patchMap_[nNewPatches] = patchi;
850  nNewPatches++;
851  }
852  }
853 
854  if (wantedPatchID == -1)
855  {
856  // Newly created patch so is at end. Check if any faces in it.
857  if (boundaryPatchSizes[oldInternalPatchID] > 0)
858  {
859  newBoundary[nNewPatches] = new internalPolyPatch
860  (
861  "oldInternalFaces",
862  boundaryPatchSizes[oldInternalPatchID],
863  patchStart,
864  nNewPatches,
865  fvMeshSubsetPtr_().boundaryMesh(),
866  internalPolyPatch::typeName
867  );
868 
869  // The index for the first patch is -1 as it originates from
870  // the internal faces
871  patchMap_[nNewPatches] = -1;
872  nNewPatches++;
873  }
874  }
875  // else
876  // {
877  // patchMap_[wantedPatchID] = -1;
878  // }
879 
880  // Reset the patch lists
881  newBoundary.setSize(nNewPatches);
882  patchMap_.setSize(nNewPatches);
883 
884  // Add the fvPatches
885  fvMeshSubsetPtr_().addFvPatches(newBoundary);
886 
887  // Subset and add any zones
888  subsetZones();
889 }
890 
891 
893 (
894  const labelList& region,
895  const label currentRegion,
896  const label patchID,
897  const bool syncPar
898 )
899 {
900  const cellList& oldCells = baseMesh().cells();
901  const faceList& oldFaces = baseMesh().faces();
902  const pointField& oldPoints = baseMesh().points();
903  const labelList& oldOwner = baseMesh().faceOwner();
904  const labelList& oldNeighbour = baseMesh().faceNeighbour();
905  const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
906  const label oldNInternalFaces = baseMesh().nInternalFaces();
907 
908  // Initial checks
909 
910  if (region.size() != oldCells.size())
911  {
913  << "Size of region " << region.size()
914  << " is not equal to number of cells in mesh " << oldCells.size()
915  << abort(FatalError);
916  }
917 
918 
919  label wantedPatchID = patchID;
920 
921  if (wantedPatchID == -1)
922  {
923  // No explicit patch specified. Put in oldInternalFaces patch.
924  // Check if patch with this name already exists.
925  wantedPatchID = oldPatches.findIndex("oldInternalFaces");
926  }
927  else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
928  {
930  << "Non-existing patch index " << wantedPatchID << endl
931  << "Should be between 0 and " << oldPatches.size()-1
932  << abort(FatalError);
933  }
934 
935  // Clear demand driven data
936  faceFlipMapPtr_.clear();
937 
938  // Get the cells for the current region.
939  cellMap_.setSize(oldCells.size());
940  label nCellsInSet = 0;
941 
942  forAll(region, oldCelli)
943  {
944  if (region[oldCelli] == currentRegion)
945  {
946  cellMap_[nCellsInSet++] = oldCelli;
947  }
948  }
949  cellMap_.setSize(nCellsInSet);
950 
951 
952  // Mark all used faces. Count number of cells using them
953  // 0: face not used anymore
954  // 1: face used by one cell, face becomes/stays boundary face
955  // 2: face still used and remains internal face
956  // 3: face coupled and used by one cell only (so should become normal,
957  // non-coupled patch face)
958  //
959  // Note that this is not really necessary - but means we can size things
960  // correctly. Also makes handling coupled faces much easier.
961 
962  labelList nCellsUsingFace(oldFaces.size(), 0);
963 
964  label nFacesInSet = 0;
965  forAll(oldFaces, oldFacei)
966  {
967  bool faceUsed = false;
968 
969  if (region[oldOwner[oldFacei]] == currentRegion)
970  {
971  nCellsUsingFace[oldFacei]++;
972  faceUsed = true;
973  }
974 
975  if
976  (
977  baseMesh().isInternalFace(oldFacei)
978  && (region[oldNeighbour[oldFacei]] == currentRegion)
979  )
980  {
981  nCellsUsingFace[oldFacei]++;
982  faceUsed = true;
983  }
984 
985  if (faceUsed)
986  {
987  nFacesInSet++;
988  }
989  }
990  faceMap_.setSize(nFacesInSet);
991 
992  // Handle coupled faces. Modifies patch faces to be uncoupled to 3.
993  Map<label> empty;
994  doCoupledPatches(syncPar, empty, nCellsUsingFace);
995 
996 
997  // See which patch to use for exposed internal faces.
998  label oldInternalPatchID = 0;
999 
1000  // Insert faces before which patch
1001  label nextPatchID = oldPatches.size();
1002 
1003  // old to new patches
1004  labelList globalPatchMap(oldPatches.size());
1005 
1006  // New patch size
1007  label nbSize = oldPatches.size();
1008 
1009  if (wantedPatchID == -1)
1010  {
1011  // Create 'oldInternalFaces' patch at the end (or before
1012  // processorPatches)
1013  // and put all exposed internal faces in there.
1014 
1015  forAll(oldPatches, patchi)
1016  {
1017  if (isA<processorPolyPatch>(oldPatches[patchi]))
1018  {
1019  nextPatchID = patchi;
1020  break;
1021  }
1022  oldInternalPatchID++;
1023  }
1024 
1025  nbSize++;
1026 
1027  // adapt old to new patches for inserted patch
1028  for (label oldPatchi = 0; oldPatchi < nextPatchID; oldPatchi++)
1029  {
1030  globalPatchMap[oldPatchi] = oldPatchi;
1031  }
1032  for
1033  (
1034  label oldPatchi = nextPatchID;
1035  oldPatchi < oldPatches.size();
1036  oldPatchi++
1037  )
1038  {
1039  globalPatchMap[oldPatchi] = oldPatchi + 1;
1040  }
1041  }
1042  else
1043  {
1044  oldInternalPatchID = wantedPatchID;
1045  nextPatchID = wantedPatchID + 1;
1046 
1047  // old to new patches
1048  globalPatchMap = identityMap(oldPatches.size());
1049  }
1050 
1051  labelList boundaryPatchSizes(nbSize, 0);
1052 
1053 
1054  // Make a global-to-local point map
1055  labelList globalPointMap(oldPoints.size(), -1);
1056 
1057  labelList globalFaceMap(oldFaces.size(), -1);
1058  label facei = 0;
1059 
1060  // 1. Pick up all preserved internal faces.
1061  for (label oldFacei = 0; oldFacei < oldNInternalFaces; oldFacei++)
1062  {
1063  if (nCellsUsingFace[oldFacei] == 2)
1064  {
1065  globalFaceMap[oldFacei] = facei;
1066  faceMap_[facei++] = oldFacei;
1067 
1068  // Mark all points from the face
1069  markPoints(oldFaces[oldFacei], globalPointMap);
1070  }
1071  }
1072 
1073  // These are all the internal faces in the mesh.
1074  label nInternalFaces = facei;
1075 
1076  // 2. Boundary faces up to where we want to insert old internal faces
1077  for
1078  (
1079  label oldPatchi = 0;
1080  oldPatchi < oldPatches.size()
1081  && oldPatchi < nextPatchID;
1082  oldPatchi++
1083  )
1084  {
1085  const polyPatch& oldPatch = oldPatches[oldPatchi];
1086 
1087  label oldFacei = oldPatch.start();
1088 
1089  forAll(oldPatch, i)
1090  {
1091  if (nCellsUsingFace[oldFacei] == 1)
1092  {
1093  // Boundary face is kept.
1094 
1095  // Mark face and increment number of points in set
1096  globalFaceMap[oldFacei] = facei;
1097  faceMap_[facei++] = oldFacei;
1098 
1099  // Mark all points from the face
1100  markPoints(oldFaces[oldFacei], globalPointMap);
1101 
1102  // Increment number of patch faces
1103  boundaryPatchSizes[globalPatchMap[oldPatchi]]++;
1104  }
1105  oldFacei++;
1106  }
1107  }
1108 
1109  // 3a. old internal faces that have become exposed.
1110  for (label oldFacei = 0; oldFacei < oldNInternalFaces; oldFacei++)
1111  {
1112  if (nCellsUsingFace[oldFacei] == 1)
1113  {
1114  globalFaceMap[oldFacei] = facei;
1115  faceMap_[facei++] = oldFacei;
1116 
1117  // Mark all points from the face
1118  markPoints(oldFaces[oldFacei], globalPointMap);
1119 
1120  // Increment number of patch faces
1121  boundaryPatchSizes[oldInternalPatchID]++;
1122  }
1123  }
1124 
1125  // 3b. coupled patch faces that have become uncoupled.
1126  for
1127  (
1128  label oldFacei = oldNInternalFaces;
1129  oldFacei < oldFaces.size();
1130  oldFacei++
1131  )
1132  {
1133  if (nCellsUsingFace[oldFacei] == 3)
1134  {
1135  globalFaceMap[oldFacei] = facei;
1136  faceMap_[facei++] = oldFacei;
1137 
1138  // Mark all points from the face
1139  markPoints(oldFaces[oldFacei], globalPointMap);
1140 
1141  // Increment number of patch faces
1142  boundaryPatchSizes[oldInternalPatchID]++;
1143  }
1144  }
1145 
1146  // 4. Remaining boundary faces
1147  for
1148  (
1149  label oldPatchi = nextPatchID;
1150  oldPatchi < oldPatches.size();
1151  oldPatchi++
1152  )
1153  {
1154  const polyPatch& oldPatch = oldPatches[oldPatchi];
1155 
1156  label oldFacei = oldPatch.start();
1157 
1158  forAll(oldPatch, i)
1159  {
1160  if (nCellsUsingFace[oldFacei] == 1)
1161  {
1162  // Boundary face is kept.
1163 
1164  // Mark face and increment number of points in set
1165  globalFaceMap[oldFacei] = facei;
1166  faceMap_[facei++] = oldFacei;
1167 
1168  // Mark all points from the face
1169  markPoints(oldFaces[oldFacei], globalPointMap);
1170 
1171  // Increment number of patch faces
1172  boundaryPatchSizes[globalPatchMap[oldPatchi]]++;
1173  }
1174  oldFacei++;
1175  }
1176  }
1177 
1178  if (facei != nFacesInSet)
1179  {
1181  << "Problem" << abort(FatalError);
1182  }
1183 
1184 
1185  // Grab the points map
1186  label nPointsInSet = 0;
1187 
1188  forAll(globalPointMap, pointi)
1189  {
1190  if (globalPointMap[pointi] != -1)
1191  {
1192  nPointsInSet++;
1193  }
1194  }
1195  pointMap_.setSize(nPointsInSet);
1196 
1197  nPointsInSet = 0;
1198 
1199  forAll(globalPointMap, pointi)
1200  {
1201  if (globalPointMap[pointi] != -1)
1202  {
1203  pointMap_[nPointsInSet] = pointi;
1204  globalPointMap[pointi] = nPointsInSet;
1205  nPointsInSet++;
1206  }
1207  }
1208 
1209 
1210  // Make a new mesh
1211  pointField newPoints(pointMap_.size());
1212 
1213  label nNewPoints = 0;
1214 
1215  forAll(pointMap_, pointi)
1216  {
1217  newPoints[nNewPoints] = oldPoints[pointMap_[pointi]];
1218  nNewPoints++;
1219  }
1220 
1221  faceList newFaces(faceMap_.size());
1222 
1223  label nNewFaces = 0;
1224 
1225  // Make internal faces
1226  for (label facei = 0; facei < nInternalFaces; facei++)
1227  {
1228  const face& oldF = oldFaces[faceMap_[facei]];
1229 
1230  face newF(oldF.size());
1231 
1232  forAll(newF, i)
1233  {
1234  newF[i] = globalPointMap[oldF[i]];
1235  }
1236 
1237  newFaces[nNewFaces] = newF;
1238  nNewFaces++;
1239  }
1240 
1241 
1242  // Make boundary faces. (different from internal since might need to be
1243  // flipped)
1244  for (label facei = nInternalFaces; facei < faceMap_.size(); facei++)
1245  {
1246  const label oldFacei = faceMap_[facei];
1247 
1248  face oldF = oldFaces[oldFacei];
1249 
1250  // Turn the faces as necessary to point outwards
1251  if (baseMesh().isInternalFace(oldFacei))
1252  {
1253  // Was internal face. Possibly turned the wrong way round
1254  if
1255  (
1256  region[oldOwner[oldFacei]] != currentRegion
1257  && region[oldNeighbour[oldFacei]] == currentRegion
1258  )
1259  {
1260  oldF = oldFaces[oldFacei].reverseFace();
1261  }
1262  }
1263 
1264  // Relabel vertices of the (possibly turned) face.
1265  face newF(oldF.size());
1266 
1267  forAll(newF, i)
1268  {
1269  newF[i] = globalPointMap[oldF[i]];
1270  }
1271 
1272  newFaces[nNewFaces] = newF;
1273  nNewFaces++;
1274  }
1275 
1276 
1277 
1278  // Create cells
1279  cellList newCells(nCellsInSet);
1280 
1281  label nNewCells = 0;
1282 
1283  forAll(cellMap_, celli)
1284  {
1285  const labelList& oldC = oldCells[cellMap_[celli]];
1286 
1287  labelList newC(oldC.size());
1288 
1289  forAll(newC, i)
1290  {
1291  newC[i] = globalFaceMap[oldC[i]];
1292  }
1293 
1294  newCells[nNewCells] = cell(newC);
1295  nNewCells++;
1296  }
1297 
1298 
1299  // Delete any old one
1300  fvMeshSubsetPtr_.clear();
1301 
1302  // Make a new mesh
1303  // Note that mesh gets registered with same name as original mesh. This is
1304  // not proper but cannot be avoided since otherwise surfaceInterpolation
1305  // cannot find its fvSchemes (it will try to read e.g.
1306  // system/region0SubSet/fvSchemes)
1307  // Make a new mesh
1308  fvMeshSubsetPtr_.reset
1309  (
1310  new fvMesh
1311  (
1312  IOobject
1313  (
1314  baseMesh().name(),
1315  baseMesh().time().name(),
1316  baseMesh().time(),
1319  ),
1320  move(newPoints),
1321  move(newFaces),
1322  move(newCells),
1323  syncPar // parallel synchronisation
1324  )
1325  );
1326 
1327  // Add old patches
1328  List<polyPatch*> newBoundary(nbSize);
1329  patchMap_.setSize(nbSize);
1330  label nNewPatches = 0;
1331  label patchStart = nInternalFaces;
1332 
1333 
1334  // For parallel: only remove patch if none of the processors has it.
1335  // This only gets done for patches before the one being inserted
1336  // (so patches < nextPatchID)
1337 
1338  // Get sum of patch sizes. Zero if patch can be deleted.
1339  labelList globalPatchSizes(boundaryPatchSizes);
1340  globalPatchSizes.setSize(nextPatchID);
1341 
1342  if (syncPar && Pstream::parRun())
1343  {
1344  // Get patch names (up to nextPatchID)
1346  patchNames[Pstream::myProcNo()] = oldPatches.names();
1347  patchNames[Pstream::myProcNo()].setSize(nextPatchID);
1350 
1351  // Get patch sizes (up to nextPatchID).
1352  // Note that up to nextPatchID the globalPatchMap is an identity so
1353  // no need to index through that.
1354  Pstream::listCombineGather(globalPatchSizes, plusEqOp<label>());
1355  Pstream::listCombineScatter(globalPatchSizes);
1356 
1357  // Now all processors have all the patchnames.
1358  // Decide: if all processors have the same patch names and size is zero
1359  // everywhere remove the patch.
1360  bool samePatches = true;
1361 
1362  for (label proci = 1; proci < patchNames.size(); proci++)
1363  {
1364  if (patchNames[proci] != patchNames[0])
1365  {
1366  samePatches = false;
1367  break;
1368  }
1369  }
1370 
1371  if (!samePatches)
1372  {
1373  // Patchnames not sync on all processors so disable removal of
1374  // zero sized patches.
1375  globalPatchSizes = labelMax;
1376  }
1377  }
1378 
1379 
1380  // Old patches
1381 
1382  for
1383  (
1384  label oldPatchi = 0;
1385  oldPatchi < oldPatches.size()
1386  && oldPatchi < nextPatchID;
1387  oldPatchi++
1388  )
1389  {
1390  const label newSize = boundaryPatchSizes[globalPatchMap[oldPatchi]];
1391 
1392  // Clone (even if 0 size)
1393  newBoundary[nNewPatches] = oldPatches[oldPatchi].clone
1394  (
1395  fvMeshSubsetPtr_().boundaryMesh(),
1396  nNewPatches,
1397  newSize,
1398  patchStart
1399  ).ptr();
1400 
1401  patchStart += newSize;
1402  patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1403  nNewPatches++;
1404  }
1405 
1406  // Inserted patch
1407 
1408  if (wantedPatchID == -1)
1409  {
1410  label oldInternalSize = boundaryPatchSizes[oldInternalPatchID];
1411 
1412  if (syncPar)
1413  {
1414  reduce(oldInternalSize, sumOp<label>());
1415  }
1416 
1417  // Newly created patch so is at end. Check if any faces in it.
1418  if (oldInternalSize > 0)
1419  {
1420  newBoundary[nNewPatches] = new internalPolyPatch
1421  (
1422  "oldInternalFaces",
1423  boundaryPatchSizes[oldInternalPatchID],
1424  patchStart,
1425  nNewPatches,
1426  fvMeshSubsetPtr_().boundaryMesh(),
1427  internalPolyPatch::typeName
1428  );
1429 
1430  // The index for the first patch is -1 as it originates from
1431  // the internal faces
1432  patchStart += boundaryPatchSizes[oldInternalPatchID];
1433  patchMap_[nNewPatches] = -1;
1434  nNewPatches++;
1435  }
1436  }
1437  // else
1438  // {
1439  // patchMap_[wantedPatchID] = -1;
1440  // }
1441 
1442  // Old patches
1443 
1444  for
1445  (
1446  label oldPatchi = nextPatchID;
1447  oldPatchi < oldPatches.size();
1448  oldPatchi++
1449  )
1450  {
1451  const label newSize = boundaryPatchSizes[globalPatchMap[oldPatchi]];
1452 
1453  // Patch still exists. Add it
1454  newBoundary[nNewPatches] = oldPatches[oldPatchi].clone
1455  (
1456  fvMeshSubsetPtr_().boundaryMesh(),
1457  nNewPatches,
1458  newSize,
1459  patchStart
1460  ).ptr();
1461 
1462  patchStart += newSize;
1463  patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1464  nNewPatches++;
1465  }
1466 
1467 
1468  // Reset the patch lists
1469  newBoundary.setSize(nNewPatches);
1470  patchMap_.setSize(nNewPatches);
1471 
1472 
1473  // Add the fvPatches
1474  fvMeshSubsetPtr_().addFvPatches(newBoundary, syncPar);
1475 
1476  // Subset and add any zones
1477  subsetZones();
1478 }
1479 
1480 
1482 (
1483  const labelHashSet& globalCellMap,
1484  const label patchID,
1485  const bool syncPar
1486 )
1487 {
1488  labelList region(baseMesh().nCells(), 0);
1489 
1490  forAllConstIter(labelHashSet, globalCellMap, iter)
1491  {
1492  region[iter.key()] = 1;
1493  }
1494  setLargeCellSubset(region, 1, patchID, syncPar);
1495 }
1496 
1497 
1499 (
1500  const labelList& region,
1501  const label currentRegion,
1502  const bool syncCouples
1503 ) const
1504 {
1505  // Collect cells to remove
1506  const labelList cellsToRemove(getCellsToRemove(region, currentRegion));
1507 
1508  return removeCells(baseMesh(), syncCouples).getExposedFaces(cellsToRemove);
1509 }
1510 
1511 
1513 (
1514  const labelList& region,
1515  const label currentRegion,
1516  const labelList& exposedFaces,
1517  const labelList& patchIDs,
1518  const bool syncCouples
1519 )
1520 {
1521  // Collect cells to remove
1522  labelList cellsToRemove(getCellsToRemove(region, currentRegion));
1523 
1524  // Mesh changing engine.
1525  polyTopoChange meshMod(baseMesh());
1526 
1527  removeCells cellRemover(baseMesh(), syncCouples);
1528 
1529  cellRemover.setRefinement
1530  (
1531  cellsToRemove,
1532  exposedFaces,
1533  patchIDs,
1534  meshMod
1535  );
1536 
1537  // Create mesh, return map from old to new mesh.
1538  autoPtr<polyTopoChangeMap> map = meshMod.makeMesh
1539  (
1540  fvMeshSubsetPtr_,
1541  IOobject
1542  (
1543  baseMesh().name(),
1544  baseMesh().time().name(),
1545  baseMesh().time(),
1548  ),
1549  baseMesh(),
1550  syncCouples
1551  );
1552 
1553  pointMap_ = map().pointMap();
1554  faceMap_ = map().faceMap();
1555  cellMap_ = map().cellMap();
1556  patchMap_ = identityMap(baseMesh().boundaryMesh().size());
1557 }
1558 
1559 
1561 {
1562  return fvMeshSubsetPtr_.valid();
1563 }
1564 
1565 
1567 {
1568  checkCellSubset();
1569 
1570  return fvMeshSubsetPtr_();
1571 }
1572 
1573 
1575 {
1576  checkCellSubset();
1577 
1578  return fvMeshSubsetPtr_();
1579 }
1580 
1581 
1583 {
1584  checkCellSubset();
1585 
1586  return pointMap_;
1587 }
1588 
1589 
1591 {
1592  checkCellSubset();
1593 
1594  return faceMap_;
1595 }
1596 
1597 
1599 {
1600  if (!faceFlipMapPtr_.valid())
1601  {
1602  const labelList& subToBaseFace = faceMap();
1603  const labelList& subToBaseCell = cellMap();
1604 
1605  faceFlipMapPtr_.reset(new labelList(subToBaseFace.size()));
1606  labelList& faceFlipMap = faceFlipMapPtr_();
1607 
1608  // Only exposed internal faces might be flipped (since we don't do
1609  // any cell renumbering, just compacting)
1610  label subInt = subMesh().nInternalFaces();
1611  const labelList& subOwn = subMesh().faceOwner();
1612  const labelList& own = baseMesh_.faceOwner();
1613 
1614  for (label subFacei = 0; subFacei < subInt; subFacei++)
1615  {
1616  faceFlipMap[subFacei] = subToBaseFace[subFacei] + 1;
1617  }
1618  for (label subFacei = subInt; subFacei < subOwn.size(); subFacei++)
1619  {
1620  const label facei = subToBaseFace[subFacei];
1621  if (subToBaseCell[subOwn[subFacei]] == own[facei])
1622  {
1623  faceFlipMap[subFacei] = facei + 1;
1624  }
1625  else
1626  {
1627  faceFlipMap[subFacei] = -facei - 1;
1628  }
1629  }
1630  }
1631 
1632  return faceFlipMapPtr_();
1633 }
1634 
1635 
1637 {
1638  checkCellSubset();
1639 
1640  return cellMap_;
1641 }
1642 
1643 
1645 {
1646  checkCellSubset();
1647 
1648  return patchMap_;
1649 }
1650 
1651 
1652 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:227
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:283
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:74
const labelList & faceMap() const
Return face map.
void setCellSubset(const labelHashSet &globalCellMap, const label patchID=-1, const bool syncPar=true)
Set the subset. Create "oldInternalFaces" patch for exposed.
Definition: fvMeshSubset.C:490
const labelList & cellMap() const
Return cell map.
fvMeshSubset(const fvMesh &)
Construct given a mesh to subset.
Definition: fvMeshSubset.C:469
bool hasSubMesh() const
Have subMesh?
const labelList & faceFlipMap() const
Return face map with sign to encode flipped faces.
virtual ~fvMeshSubset()
Destructor.
Definition: fvMeshSubset.C:483
void setLargeCellSubset(const labelList &region, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
Set the subset from all cells with region == currentRegion.
Definition: fvMeshSubset.C:893
const fvMesh & subMesh() const
Return reference to subset mesh.
const labelList & patchMap() const
Return patch map.
const labelList & pointMap() const
Return point map.
labelList getExposedFaces(const labelList &region, const label currentRegion, const bool syncCouples=true) const
Get labels of exposed faces.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Constraint patch to hold internal faces exposed by sub-setting.
Foam::polyBoundaryMesh.
label findIndex(const word &patchName) const
Find patch index given a name.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
wordList names() const
Return the list of patch names.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
Direct mesh changes based on v1.3 polyTopoChange syntax.
autoPtr< polyTopoChangeMap > makeMesh(autoPtr< fvMesh > &newMesh, const IOobject &io, const polyMesh &mesh, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Create new mesh with old mesh patches.
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:60
void setRefinement(const labelList &cellsToRemove, const labelList &facesToExpose, const labelList &patchIDs, polyTopoChange &) const
Play commands into polyTopoChange to remove cells.
Definition: removeCells.C:175
labelList getExposedFaces(const labelList &cellsToRemove) const
Get labels of exposed faces.
Definition: removeCells.C:69
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
label nPoints
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void sort(UList< T > &)
Definition: UList.C:115
defineTypeNameAndDebug(combustionModel, 0)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
error FatalError
static const label labelMax
Definition: label.H:62
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
wordList patchNames(nPatches)