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-2023 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 meshPointZones& 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  i,
341  fvMeshSubsetPtr_().pointZones()
342  );
343  }
344 
345 
346  // FaceZones
347  const meshFaceZones& faceZones = baseMesh().faceZones();
348 
349  // Do we need to remove zones where the side we're interested in
350  // no longer exists? Guess not.
351  List<faceZone*> fZonePtrs(faceZones.size());
352 
353  forAll(faceZones, i)
354  {
355  const faceZone& fz = faceZones[i];
356 
357  // Expand faceZone to full mesh
358  // +1 : part of faceZone, flipped
359  // -1 : ,, , unflipped
360  // 0 : not part of faceZone
361  labelList zone(baseMesh().nFaces(), 0);
362  forAll(fz, j)
363  {
364  if (fz.flipMap()[j])
365  {
366  zone[fz[j]] = 1;
367  }
368  else
369  {
370  zone[fz[j]] = -1;
371  }
372  }
373 
374  // Select faces
375  label nSub = 0;
376  forAll(faceMap(), j)
377  {
378  if (zone[faceMap()[j]] != 0)
379  {
380  nSub++;
381  }
382  }
383  labelList subAddressing(nSub);
384  boolList subFlipStatus(nSub);
385  nSub = 0;
386  forAll(faceMap(), subFacei)
387  {
388  const label meshFacei = faceMap()[subFacei];
389  if (zone[meshFacei] != 0)
390  {
391  subAddressing[nSub] = subFacei;
392  const label subOwner = subMesh().faceOwner()[subFacei];
393  const label baseOwner = baseMesh().faceOwner()[meshFacei];
394  // If subowner is the same cell as the base keep the flip status
395  const bool sameOwner = (cellMap()[subOwner] == baseOwner);
396  const bool flip = (zone[meshFacei] == 1);
397  subFlipStatus[nSub] = (sameOwner == flip);
398 
399  nSub++;
400  }
401  }
402 
403  fZonePtrs[i] = new faceZone
404  (
405  fz.name(),
406  subAddressing,
407  subFlipStatus,
408  i,
409  fvMeshSubsetPtr_().faceZones()
410  );
411  }
412 
413  // CellZones
414  const meshCellZones& cellZones = baseMesh().cellZones();
415 
416  List<cellZone*> cZonePtrs(cellZones.size());
417 
418  forAll(cellZones, i)
419  {
420  const cellZone& cz = cellZones[i];
421 
422  cZonePtrs[i] = new cellZone
423  (
424  cz.name(),
425  subset(baseMesh().nCells(), cz, cellMap()),
426  i,
427  fvMeshSubsetPtr_().cellZones()
428  );
429  }
430 
431 
432  // Add the zones
433  fvMeshSubsetPtr_().addZones(pZonePtrs, fZonePtrs, cZonePtrs);
434 }
435 
436 
437 Foam::labelList Foam::fvMeshSubset::getCellsToRemove
438 (
439  const labelList& region,
440  const label currentRegion
441 ) const
442 {
443  // Count
444  label nKeep = 0;
445  forAll(region, cellI)
446  {
447  if (region[cellI] == currentRegion)
448  {
449  nKeep++;
450  }
451  }
452 
453  // Collect cells to remove
454  label nRemove = baseMesh().nCells() - nKeep;
455  labelList cellsToRemove(nRemove);
456 
457  nRemove = 0;
458  forAll(region, cellI)
459  {
460  if (region[cellI] != currentRegion)
461  {
462  cellsToRemove[nRemove++] = cellI;
463  }
464  }
465 
466  return cellsToRemove;
467 }
468 
469 
470 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
471 
473 :
474  baseMesh_(baseMesh),
475  fvMeshSubsetPtr_(nullptr),
476  pointMap_(0),
477  faceMap_(0),
478  cellMap_(0),
479  patchMap_(0),
480  faceFlipMapPtr_()
481 {}
482 
483 
484 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
485 
487 {}
488 
489 
490 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
491 
493 (
494  const labelHashSet& globalCellMap,
495  const label patchID,
496  const bool syncPar
497 )
498 {
499  // Initial check on patches before doing anything time consuming.
500  const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
501  const cellList& oldCells = baseMesh().cells();
502  const faceList& oldFaces = baseMesh().faces();
503  const pointField& oldPoints = baseMesh().points();
504  const labelList& oldOwner = baseMesh().faceOwner();
505  const labelList& oldNeighbour = baseMesh().faceNeighbour();
506 
507  label wantedPatchID = patchID;
508 
509  if (wantedPatchID == -1)
510  {
511  // No explicit patch specified. Put in oldInternalFaces patch.
512  // Check if patch with this name already exists.
513  wantedPatchID = oldPatches.findPatchID("oldInternalFaces");
514  }
515  else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
516  {
518  << "Non-existing patch index " << wantedPatchID << endl
519  << "Should be between 0 and " << oldPatches.size()-1
520  << abort(FatalError);
521  }
522 
523 
524  // Clear demand driven data
525  faceFlipMapPtr_.clear();
526 
527 
528  cellMap_ = globalCellMap.toc();
529 
530  // Sort the cell map in the ascending order
531  sort(cellMap_);
532 
533  // Approximate sizing parameters for face and point lists
534  const label avgNFacesPerCell = 6;
535  const label avgNPointsPerFace = 4;
536 
537 
538  const label nCellsInSet = cellMap_.size();
539 
540  // Mark all used faces
541 
542  Map<label> facesToSubset(avgNFacesPerCell*nCellsInSet);
543 
544  forAll(cellMap_, celli)
545  {
546  // Mark all faces from the cell
547  const labelList& curFaces = oldCells[cellMap_[celli]];
548 
549  forAll(curFaces, facei)
550  {
551  if (!facesToSubset.found(curFaces[facei]))
552  {
553  facesToSubset.insert(curFaces[facei], 1);
554  }
555  else
556  {
557  facesToSubset[curFaces[facei]]++;
558  }
559  }
560  }
561 
562  // Handle coupled faces. Modifies patch faces to be uncoupled to 3.
563  labelList empty;
564  doCoupledPatches(syncPar, facesToSubset, empty);
565 
566  // Mark all used points and make a global-to-local face map
567  Map<label> globalFaceMap(facesToSubset.size());
568 
569  // Make a global-to-local point map
570  Map<label> globalPointMap(avgNPointsPerFace*facesToSubset.size());
571 
572  // This is done in two goes, so that the boundary faces are last
573  // in the list. Because of this, I need to create the face map
574  // along the way rather than just grab the table of contents.
575  labelList facesToc = facesToSubset.toc();
576  sort(facesToc);
577  faceMap_.setSize(facesToc.size());
578 
579  // 1. Get all faces that will be internal to the submesh.
580  forAll(facesToc, facei)
581  {
582  if (facesToSubset[facesToc[facei]] == 2)
583  {
584  // Mark face and increment number of points in set
585  faceMap_[globalFaceMap.size()] = facesToc[facei];
586  globalFaceMap.insert(facesToc[facei], globalFaceMap.size());
587 
588  // Mark all points from the face
589  markPoints(oldFaces[facesToc[facei]], globalPointMap);
590  }
591  }
592 
593  // These are all the internal faces in the mesh.
594  const label nInternalFaces = globalFaceMap.size();
595 
596 
597  // Where to insert old internal faces.
598  label oldPatchStart = labelMax;
599  if (wantedPatchID != -1)
600  {
601  oldPatchStart = oldPatches[wantedPatchID].start();
602  }
603 
604 
605  label facei = 0;
606 
607  // 2. Boundary faces up to where we want to insert old internal faces
608  for (; facei< facesToc.size(); facei++)
609  {
610  if (facesToc[facei] >= oldPatchStart)
611  {
612  break;
613  }
614  if
615  (
616  !baseMesh().isInternalFace(facesToc[facei])
617  && facesToSubset[facesToc[facei]] == 1
618  )
619  {
620  // Mark face and increment number of points in set
621  faceMap_[globalFaceMap.size()] = facesToc[facei];
622  globalFaceMap.insert(facesToc[facei], globalFaceMap.size());
623 
624  // Mark all points from the face
625  markPoints(oldFaces[facesToc[facei]], globalPointMap);
626  }
627  }
628 
629  // 3. old internal faces and uncoupled faces
630  forAll(facesToc, intFacei)
631  {
632  if
633  (
634  (
635  baseMesh().isInternalFace(facesToc[intFacei])
636  && facesToSubset[facesToc[intFacei]] == 1
637  )
638  || (
639  !baseMesh().isInternalFace(facesToc[intFacei])
640  && facesToSubset[facesToc[intFacei]] == 3
641  )
642  )
643  {
644  // Mark face and increment number of points in set
645  faceMap_[globalFaceMap.size()] = facesToc[intFacei];
646  globalFaceMap.insert(facesToc[intFacei], globalFaceMap.size());
647 
648  // Mark all points from the face
649  markPoints(oldFaces[facesToc[intFacei]], globalPointMap);
650  }
651  }
652 
653  // 4. Remaining boundary faces
654  for (; facei< facesToc.size(); facei++)
655  {
656  if
657  (
658  !baseMesh().isInternalFace(facesToc[facei])
659  && facesToSubset[facesToc[facei]] == 1
660  )
661  {
662  // Mark face and increment number of points in set
663  faceMap_[globalFaceMap.size()] = facesToc[facei];
664  globalFaceMap.insert(facesToc[facei], globalFaceMap.size());
665 
666  // Mark all points from the face
667  markPoints(oldFaces[facesToc[facei]], globalPointMap);
668  }
669  }
670 
671 
672 
673  // Grab the points map
674  pointMap_ = globalPointMap.toc();
675  sort(pointMap_);
676 
677  forAll(pointMap_, pointi)
678  {
679  globalPointMap[pointMap_[pointi]] = pointi;
680  }
681 
682  // Make a new mesh
683  pointField newPoints(globalPointMap.size());
684 
685  label nNewPoints = 0;
686 
687  forAll(pointMap_, pointi)
688  {
689  newPoints[nNewPoints] = oldPoints[pointMap_[pointi]];
690  nNewPoints++;
691  }
692 
693  faceList newFaces(globalFaceMap.size());
694 
695  label nNewFaces = 0;
696 
697  // Make internal faces
698  for (label facei = 0; facei < nInternalFaces; facei++)
699  {
700  const face& oldF = oldFaces[faceMap_[facei]];
701 
702  face newF(oldF.size());
703 
704  forAll(newF, i)
705  {
706  newF[i] = globalPointMap[oldF[i]];
707  }
708 
709  newFaces[nNewFaces] = newF;
710  nNewFaces++;
711  }
712 
713  // Make boundary faces
714 
715  label nbSize = oldPatches.size();
716  label oldInternalPatchID = -1;
717 
718  if (wantedPatchID == -1)
719  {
720  // Create 'oldInternalFaces' patch at the end
721  // and put all exposed internal faces in there.
722  oldInternalPatchID = nbSize;
723  nbSize++;
724 
725  }
726  else
727  {
728  oldInternalPatchID = wantedPatchID;
729  }
730 
731 
732  // Grad size and start of each patch on the fly. Because of the
733  // structure of the underlying mesh, the patches will appear in the
734  // ascending order
735  labelList boundaryPatchSizes(nbSize, 0);
736 
737  // Assign boundary faces. Visited in order of faceMap_.
738  for (label facei = nInternalFaces; facei < faceMap_.size(); facei++)
739  {
740  const label oldFacei = faceMap_[facei];
741 
742  face oldF = oldFaces[oldFacei];
743 
744  // Turn the faces as necessary to point outwards
745  if (baseMesh().isInternalFace(oldFacei))
746  {
747  // Internal face. Possibly turned the wrong way round
748  if
749  (
750  !globalCellMap.found(oldOwner[oldFacei])
751  && globalCellMap.found(oldNeighbour[oldFacei])
752  )
753  {
754  oldF = oldFaces[oldFacei].reverseFace();
755  }
756 
757  // Update count for patch
758  boundaryPatchSizes[oldInternalPatchID]++;
759  }
760  else if (facesToSubset[oldFacei] == 3)
761  {
762  // Uncoupled face. Increment the old patch.
763  boundaryPatchSizes[oldInternalPatchID]++;
764  }
765  else
766  {
767  // Boundary face. Increment the appropriate patch
768  const label patchOfFace = oldPatches.whichPatch(oldFacei);
769 
770  // Update count for patch
771  boundaryPatchSizes[patchOfFace]++;
772  }
773 
774  face newF(oldF.size());
775 
776  forAll(newF, i)
777  {
778  newF[i] = globalPointMap[oldF[i]];
779  }
780 
781  newFaces[nNewFaces] = newF;
782  nNewFaces++;
783  }
784 
785 
786 
787  // Create cells
788  cellList newCells(nCellsInSet);
789 
790  label nNewCells = 0;
791 
792  forAll(cellMap_, celli)
793  {
794  const labelList& oldC = oldCells[cellMap_[celli]];
795 
796  labelList newC(oldC.size());
797 
798  forAll(newC, i)
799  {
800  newC[i] = globalFaceMap[oldC[i]];
801  }
802 
803  newCells[nNewCells] = cell(newC);
804  nNewCells++;
805  }
806 
807 
808  // Delete any old mesh
809  fvMeshSubsetPtr_.clear();
810 
811  // Make a new mesh
812  fvMeshSubsetPtr_.reset
813  (
814  new fvMesh
815  (
816  IOobject
817  (
818  baseMesh().name() + "SubSet",
819  baseMesh().time().name(),
820  baseMesh().time(),
823  ),
824  move(newPoints),
825  move(newFaces),
826  move(newCells)
827  )
828  );
829 
830 
831  // Add old patches
832  List<polyPatch*> newBoundary(nbSize);
833  patchMap_.setSize(nbSize);
834  label nNewPatches = 0;
835  label patchStart = nInternalFaces;
836 
837 
838  forAll(oldPatches, patchi)
839  {
840  if (boundaryPatchSizes[patchi] > 0)
841  {
842  // Patch still exists. Add it
843  newBoundary[nNewPatches] = oldPatches[patchi].clone
844  (
845  fvMeshSubsetPtr_().boundaryMesh(),
846  nNewPatches,
847  boundaryPatchSizes[patchi],
848  patchStart
849  ).ptr();
850 
851  patchStart += boundaryPatchSizes[patchi];
852  patchMap_[nNewPatches] = patchi;
853  nNewPatches++;
854  }
855  }
856 
857  if (wantedPatchID == -1)
858  {
859  // Newly created patch so is at end. Check if any faces in it.
860  if (boundaryPatchSizes[oldInternalPatchID] > 0)
861  {
862  newBoundary[nNewPatches] = new internalPolyPatch
863  (
864  "oldInternalFaces",
865  boundaryPatchSizes[oldInternalPatchID],
866  patchStart,
867  nNewPatches,
868  fvMeshSubsetPtr_().boundaryMesh(),
869  internalPolyPatch::typeName
870  );
871 
872  // The index for the first patch is -1 as it originates from
873  // the internal faces
874  patchMap_[nNewPatches] = -1;
875  nNewPatches++;
876  }
877  }
878  // else
879  // {
880  // patchMap_[wantedPatchID] = -1;
881  // }
882 
883  // Reset the patch lists
884  newBoundary.setSize(nNewPatches);
885  patchMap_.setSize(nNewPatches);
886 
887  // Add the fvPatches
888  fvMeshSubsetPtr_().addFvPatches(newBoundary);
889 
890  // Subset and add any zones
891  subsetZones();
892 }
893 
894 
896 (
897  const labelList& region,
898  const label currentRegion,
899  const label patchID,
900  const bool syncPar
901 )
902 {
903  const cellList& oldCells = baseMesh().cells();
904  const faceList& oldFaces = baseMesh().faces();
905  const pointField& oldPoints = baseMesh().points();
906  const labelList& oldOwner = baseMesh().faceOwner();
907  const labelList& oldNeighbour = baseMesh().faceNeighbour();
908  const polyBoundaryMesh& oldPatches = baseMesh().boundaryMesh();
909  const label oldNInternalFaces = baseMesh().nInternalFaces();
910 
911  // Initial checks
912 
913  if (region.size() != oldCells.size())
914  {
916  << "Size of region " << region.size()
917  << " is not equal to number of cells in mesh " << oldCells.size()
918  << abort(FatalError);
919  }
920 
921 
922  label wantedPatchID = patchID;
923 
924  if (wantedPatchID == -1)
925  {
926  // No explicit patch specified. Put in oldInternalFaces patch.
927  // Check if patch with this name already exists.
928  wantedPatchID = oldPatches.findPatchID("oldInternalFaces");
929  }
930  else if (wantedPatchID < 0 || wantedPatchID >= oldPatches.size())
931  {
933  << "Non-existing patch index " << wantedPatchID << endl
934  << "Should be between 0 and " << oldPatches.size()-1
935  << abort(FatalError);
936  }
937 
938  // Clear demand driven data
939  faceFlipMapPtr_.clear();
940 
941  // Get the cells for the current region.
942  cellMap_.setSize(oldCells.size());
943  label nCellsInSet = 0;
944 
945  forAll(region, oldCelli)
946  {
947  if (region[oldCelli] == currentRegion)
948  {
949  cellMap_[nCellsInSet++] = oldCelli;
950  }
951  }
952  cellMap_.setSize(nCellsInSet);
953 
954 
955  // Mark all used faces. Count number of cells using them
956  // 0: face not used anymore
957  // 1: face used by one cell, face becomes/stays boundary face
958  // 2: face still used and remains internal face
959  // 3: face coupled and used by one cell only (so should become normal,
960  // non-coupled patch face)
961  //
962  // Note that this is not really necessary - but means we can size things
963  // correctly. Also makes handling coupled faces much easier.
964 
965  labelList nCellsUsingFace(oldFaces.size(), 0);
966 
967  label nFacesInSet = 0;
968  forAll(oldFaces, oldFacei)
969  {
970  bool faceUsed = false;
971 
972  if (region[oldOwner[oldFacei]] == currentRegion)
973  {
974  nCellsUsingFace[oldFacei]++;
975  faceUsed = true;
976  }
977 
978  if
979  (
980  baseMesh().isInternalFace(oldFacei)
981  && (region[oldNeighbour[oldFacei]] == currentRegion)
982  )
983  {
984  nCellsUsingFace[oldFacei]++;
985  faceUsed = true;
986  }
987 
988  if (faceUsed)
989  {
990  nFacesInSet++;
991  }
992  }
993  faceMap_.setSize(nFacesInSet);
994 
995  // Handle coupled faces. Modifies patch faces to be uncoupled to 3.
996  Map<label> empty;
997  doCoupledPatches(syncPar, empty, nCellsUsingFace);
998 
999 
1000  // See which patch to use for exposed internal faces.
1001  label oldInternalPatchID = 0;
1002 
1003  // Insert faces before which patch
1004  label nextPatchID = oldPatches.size();
1005 
1006  // old to new patches
1007  labelList globalPatchMap(oldPatches.size());
1008 
1009  // New patch size
1010  label nbSize = oldPatches.size();
1011 
1012  if (wantedPatchID == -1)
1013  {
1014  // Create 'oldInternalFaces' patch at the end (or before
1015  // processorPatches)
1016  // and put all exposed internal faces in there.
1017 
1018  forAll(oldPatches, patchi)
1019  {
1020  if (isA<processorPolyPatch>(oldPatches[patchi]))
1021  {
1022  nextPatchID = patchi;
1023  break;
1024  }
1025  oldInternalPatchID++;
1026  }
1027 
1028  nbSize++;
1029 
1030  // adapt old to new patches for inserted patch
1031  for (label oldPatchi = 0; oldPatchi < nextPatchID; oldPatchi++)
1032  {
1033  globalPatchMap[oldPatchi] = oldPatchi;
1034  }
1035  for
1036  (
1037  label oldPatchi = nextPatchID;
1038  oldPatchi < oldPatches.size();
1039  oldPatchi++
1040  )
1041  {
1042  globalPatchMap[oldPatchi] = oldPatchi + 1;
1043  }
1044  }
1045  else
1046  {
1047  oldInternalPatchID = wantedPatchID;
1048  nextPatchID = wantedPatchID + 1;
1049 
1050  // old to new patches
1051  globalPatchMap = identityMap(oldPatches.size());
1052  }
1053 
1054  labelList boundaryPatchSizes(nbSize, 0);
1055 
1056 
1057  // Make a global-to-local point map
1058  labelList globalPointMap(oldPoints.size(), -1);
1059 
1060  labelList globalFaceMap(oldFaces.size(), -1);
1061  label facei = 0;
1062 
1063  // 1. Pick up all preserved internal faces.
1064  for (label oldFacei = 0; oldFacei < oldNInternalFaces; oldFacei++)
1065  {
1066  if (nCellsUsingFace[oldFacei] == 2)
1067  {
1068  globalFaceMap[oldFacei] = facei;
1069  faceMap_[facei++] = oldFacei;
1070 
1071  // Mark all points from the face
1072  markPoints(oldFaces[oldFacei], globalPointMap);
1073  }
1074  }
1075 
1076  // These are all the internal faces in the mesh.
1077  label nInternalFaces = facei;
1078 
1079  // 2. Boundary faces up to where we want to insert old internal faces
1080  for
1081  (
1082  label oldPatchi = 0;
1083  oldPatchi < oldPatches.size()
1084  && oldPatchi < nextPatchID;
1085  oldPatchi++
1086  )
1087  {
1088  const polyPatch& oldPatch = oldPatches[oldPatchi];
1089 
1090  label oldFacei = oldPatch.start();
1091 
1092  forAll(oldPatch, i)
1093  {
1094  if (nCellsUsingFace[oldFacei] == 1)
1095  {
1096  // Boundary face is kept.
1097 
1098  // Mark face and increment number of points in set
1099  globalFaceMap[oldFacei] = facei;
1100  faceMap_[facei++] = oldFacei;
1101 
1102  // Mark all points from the face
1103  markPoints(oldFaces[oldFacei], globalPointMap);
1104 
1105  // Increment number of patch faces
1106  boundaryPatchSizes[globalPatchMap[oldPatchi]]++;
1107  }
1108  oldFacei++;
1109  }
1110  }
1111 
1112  // 3a. old internal faces that have become exposed.
1113  for (label oldFacei = 0; oldFacei < oldNInternalFaces; oldFacei++)
1114  {
1115  if (nCellsUsingFace[oldFacei] == 1)
1116  {
1117  globalFaceMap[oldFacei] = facei;
1118  faceMap_[facei++] = oldFacei;
1119 
1120  // Mark all points from the face
1121  markPoints(oldFaces[oldFacei], globalPointMap);
1122 
1123  // Increment number of patch faces
1124  boundaryPatchSizes[oldInternalPatchID]++;
1125  }
1126  }
1127 
1128  // 3b. coupled patch faces that have become uncoupled.
1129  for
1130  (
1131  label oldFacei = oldNInternalFaces;
1132  oldFacei < oldFaces.size();
1133  oldFacei++
1134  )
1135  {
1136  if (nCellsUsingFace[oldFacei] == 3)
1137  {
1138  globalFaceMap[oldFacei] = facei;
1139  faceMap_[facei++] = oldFacei;
1140 
1141  // Mark all points from the face
1142  markPoints(oldFaces[oldFacei], globalPointMap);
1143 
1144  // Increment number of patch faces
1145  boundaryPatchSizes[oldInternalPatchID]++;
1146  }
1147  }
1148 
1149  // 4. Remaining boundary faces
1150  for
1151  (
1152  label oldPatchi = nextPatchID;
1153  oldPatchi < oldPatches.size();
1154  oldPatchi++
1155  )
1156  {
1157  const polyPatch& oldPatch = oldPatches[oldPatchi];
1158 
1159  label oldFacei = oldPatch.start();
1160 
1161  forAll(oldPatch, i)
1162  {
1163  if (nCellsUsingFace[oldFacei] == 1)
1164  {
1165  // Boundary face is kept.
1166 
1167  // Mark face and increment number of points in set
1168  globalFaceMap[oldFacei] = facei;
1169  faceMap_[facei++] = oldFacei;
1170 
1171  // Mark all points from the face
1172  markPoints(oldFaces[oldFacei], globalPointMap);
1173 
1174  // Increment number of patch faces
1175  boundaryPatchSizes[globalPatchMap[oldPatchi]]++;
1176  }
1177  oldFacei++;
1178  }
1179  }
1180 
1181  if (facei != nFacesInSet)
1182  {
1184  << "Problem" << abort(FatalError);
1185  }
1186 
1187 
1188  // Grab the points map
1189  label nPointsInSet = 0;
1190 
1191  forAll(globalPointMap, pointi)
1192  {
1193  if (globalPointMap[pointi] != -1)
1194  {
1195  nPointsInSet++;
1196  }
1197  }
1198  pointMap_.setSize(nPointsInSet);
1199 
1200  nPointsInSet = 0;
1201 
1202  forAll(globalPointMap, pointi)
1203  {
1204  if (globalPointMap[pointi] != -1)
1205  {
1206  pointMap_[nPointsInSet] = pointi;
1207  globalPointMap[pointi] = nPointsInSet;
1208  nPointsInSet++;
1209  }
1210  }
1211 
1212 
1213  // Make a new mesh
1214  pointField newPoints(pointMap_.size());
1215 
1216  label nNewPoints = 0;
1217 
1218  forAll(pointMap_, pointi)
1219  {
1220  newPoints[nNewPoints] = oldPoints[pointMap_[pointi]];
1221  nNewPoints++;
1222  }
1223 
1224  faceList newFaces(faceMap_.size());
1225 
1226  label nNewFaces = 0;
1227 
1228  // Make internal faces
1229  for (label facei = 0; facei < nInternalFaces; facei++)
1230  {
1231  const face& oldF = oldFaces[faceMap_[facei]];
1232 
1233  face newF(oldF.size());
1234 
1235  forAll(newF, i)
1236  {
1237  newF[i] = globalPointMap[oldF[i]];
1238  }
1239 
1240  newFaces[nNewFaces] = newF;
1241  nNewFaces++;
1242  }
1243 
1244 
1245  // Make boundary faces. (different from internal since might need to be
1246  // flipped)
1247  for (label facei = nInternalFaces; facei < faceMap_.size(); facei++)
1248  {
1249  const label oldFacei = faceMap_[facei];
1250 
1251  face oldF = oldFaces[oldFacei];
1252 
1253  // Turn the faces as necessary to point outwards
1254  if (baseMesh().isInternalFace(oldFacei))
1255  {
1256  // Was internal face. Possibly turned the wrong way round
1257  if
1258  (
1259  region[oldOwner[oldFacei]] != currentRegion
1260  && region[oldNeighbour[oldFacei]] == currentRegion
1261  )
1262  {
1263  oldF = oldFaces[oldFacei].reverseFace();
1264  }
1265  }
1266 
1267  // Relabel vertices of the (possibly turned) face.
1268  face newF(oldF.size());
1269 
1270  forAll(newF, i)
1271  {
1272  newF[i] = globalPointMap[oldF[i]];
1273  }
1274 
1275  newFaces[nNewFaces] = newF;
1276  nNewFaces++;
1277  }
1278 
1279 
1280 
1281  // Create cells
1282  cellList newCells(nCellsInSet);
1283 
1284  label nNewCells = 0;
1285 
1286  forAll(cellMap_, celli)
1287  {
1288  const labelList& oldC = oldCells[cellMap_[celli]];
1289 
1290  labelList newC(oldC.size());
1291 
1292  forAll(newC, i)
1293  {
1294  newC[i] = globalFaceMap[oldC[i]];
1295  }
1296 
1297  newCells[nNewCells] = cell(newC);
1298  nNewCells++;
1299  }
1300 
1301 
1302  // Delete any old one
1303  fvMeshSubsetPtr_.clear();
1304 
1305  // Make a new mesh
1306  // Note that mesh gets registered with same name as original mesh. This is
1307  // not proper but cannot be avoided since otherwise surfaceInterpolation
1308  // cannot find its fvSchemes (it will try to read e.g.
1309  // system/region0SubSet/fvSchemes)
1310  // Make a new mesh
1311  fvMeshSubsetPtr_.reset
1312  (
1313  new fvMesh
1314  (
1315  IOobject
1316  (
1317  baseMesh().name(),
1318  baseMesh().time().name(),
1319  baseMesh().time(),
1322  ),
1323  move(newPoints),
1324  move(newFaces),
1325  move(newCells),
1326  syncPar // parallel synchronisation
1327  )
1328  );
1329 
1330  // Add old patches
1331  List<polyPatch*> newBoundary(nbSize);
1332  patchMap_.setSize(nbSize);
1333  label nNewPatches = 0;
1334  label patchStart = nInternalFaces;
1335 
1336 
1337  // For parallel: only remove patch if none of the processors has it.
1338  // This only gets done for patches before the one being inserted
1339  // (so patches < nextPatchID)
1340 
1341  // Get sum of patch sizes. Zero if patch can be deleted.
1342  labelList globalPatchSizes(boundaryPatchSizes);
1343  globalPatchSizes.setSize(nextPatchID);
1344 
1345  if (syncPar && Pstream::parRun())
1346  {
1347  // Get patch names (up to nextPatchID)
1349  patchNames[Pstream::myProcNo()] = oldPatches.names();
1350  patchNames[Pstream::myProcNo()].setSize(nextPatchID);
1353 
1354  // Get patch sizes (up to nextPatchID).
1355  // Note that up to nextPatchID the globalPatchMap is an identity so
1356  // no need to index through that.
1357  Pstream::listCombineGather(globalPatchSizes, plusEqOp<label>());
1358  Pstream::listCombineScatter(globalPatchSizes);
1359 
1360  // Now all processors have all the patchnames.
1361  // Decide: if all processors have the same patch names and size is zero
1362  // everywhere remove the patch.
1363  bool samePatches = true;
1364 
1365  for (label proci = 1; proci < patchNames.size(); proci++)
1366  {
1367  if (patchNames[proci] != patchNames[0])
1368  {
1369  samePatches = false;
1370  break;
1371  }
1372  }
1373 
1374  if (!samePatches)
1375  {
1376  // Patchnames not sync on all processors so disable removal of
1377  // zero sized patches.
1378  globalPatchSizes = labelMax;
1379  }
1380  }
1381 
1382 
1383  // Old patches
1384 
1385  for
1386  (
1387  label oldPatchi = 0;
1388  oldPatchi < oldPatches.size()
1389  && oldPatchi < nextPatchID;
1390  oldPatchi++
1391  )
1392  {
1393  const label newSize = boundaryPatchSizes[globalPatchMap[oldPatchi]];
1394 
1395  // Clone (even if 0 size)
1396  newBoundary[nNewPatches] = oldPatches[oldPatchi].clone
1397  (
1398  fvMeshSubsetPtr_().boundaryMesh(),
1399  nNewPatches,
1400  newSize,
1401  patchStart
1402  ).ptr();
1403 
1404  patchStart += newSize;
1405  patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1406  nNewPatches++;
1407  }
1408 
1409  // Inserted patch
1410 
1411  if (wantedPatchID == -1)
1412  {
1413  label oldInternalSize = boundaryPatchSizes[oldInternalPatchID];
1414 
1415  if (syncPar)
1416  {
1417  reduce(oldInternalSize, sumOp<label>());
1418  }
1419 
1420  // Newly created patch so is at end. Check if any faces in it.
1421  if (oldInternalSize > 0)
1422  {
1423  newBoundary[nNewPatches] = new internalPolyPatch
1424  (
1425  "oldInternalFaces",
1426  boundaryPatchSizes[oldInternalPatchID],
1427  patchStart,
1428  nNewPatches,
1429  fvMeshSubsetPtr_().boundaryMesh(),
1430  internalPolyPatch::typeName
1431  );
1432 
1433  // The index for the first patch is -1 as it originates from
1434  // the internal faces
1435  patchStart += boundaryPatchSizes[oldInternalPatchID];
1436  patchMap_[nNewPatches] = -1;
1437  nNewPatches++;
1438  }
1439  }
1440  // else
1441  // {
1442  // patchMap_[wantedPatchID] = -1;
1443  // }
1444 
1445  // Old patches
1446 
1447  for
1448  (
1449  label oldPatchi = nextPatchID;
1450  oldPatchi < oldPatches.size();
1451  oldPatchi++
1452  )
1453  {
1454  const label newSize = boundaryPatchSizes[globalPatchMap[oldPatchi]];
1455 
1456  // Patch still exists. Add it
1457  newBoundary[nNewPatches] = oldPatches[oldPatchi].clone
1458  (
1459  fvMeshSubsetPtr_().boundaryMesh(),
1460  nNewPatches,
1461  newSize,
1462  patchStart
1463  ).ptr();
1464 
1465  patchStart += newSize;
1466  patchMap_[nNewPatches] = oldPatchi; // compact patchMap
1467  nNewPatches++;
1468  }
1469 
1470 
1471  // Reset the patch lists
1472  newBoundary.setSize(nNewPatches);
1473  patchMap_.setSize(nNewPatches);
1474 
1475 
1476  // Add the fvPatches
1477  fvMeshSubsetPtr_().addFvPatches(newBoundary, syncPar);
1478 
1479  // Subset and add any zones
1480  subsetZones();
1481 }
1482 
1483 
1485 (
1486  const labelHashSet& globalCellMap,
1487  const label patchID,
1488  const bool syncPar
1489 )
1490 {
1491  labelList region(baseMesh().nCells(), 0);
1492 
1493  forAllConstIter(labelHashSet, globalCellMap, iter)
1494  {
1495  region[iter.key()] = 1;
1496  }
1497  setLargeCellSubset(region, 1, patchID, syncPar);
1498 }
1499 
1500 
1502 (
1503  const labelList& region,
1504  const label currentRegion,
1505  const bool syncCouples
1506 ) const
1507 {
1508  // Collect cells to remove
1509  const labelList cellsToRemove(getCellsToRemove(region, currentRegion));
1510 
1511  return removeCells(baseMesh(), syncCouples).getExposedFaces(cellsToRemove);
1512 }
1513 
1514 
1516 (
1517  const labelList& region,
1518  const label currentRegion,
1519  const labelList& exposedFaces,
1520  const labelList& patchIDs,
1521  const bool syncCouples
1522 )
1523 {
1524  // Collect cells to remove
1525  labelList cellsToRemove(getCellsToRemove(region, currentRegion));
1526 
1527  // Mesh changing engine.
1528  polyTopoChange meshMod(baseMesh());
1529 
1530  removeCells cellRemover(baseMesh(), syncCouples);
1531 
1532  cellRemover.setRefinement
1533  (
1534  cellsToRemove,
1535  exposedFaces,
1536  patchIDs,
1537  meshMod
1538  );
1539 
1540  // Create mesh, return map from old to new mesh.
1541  autoPtr<polyTopoChangeMap> map = meshMod.makeMesh
1542  (
1543  fvMeshSubsetPtr_,
1544  IOobject
1545  (
1546  baseMesh().name(),
1547  baseMesh().time().name(),
1548  baseMesh().time(),
1551  ),
1552  baseMesh(),
1553  syncCouples
1554  );
1555 
1556  pointMap_ = map().pointMap();
1557  faceMap_ = map().faceMap();
1558  cellMap_ = map().cellMap();
1559  patchMap_ = identityMap(baseMesh().boundaryMesh().size());
1560 }
1561 
1562 
1564 {
1565  return fvMeshSubsetPtr_.valid();
1566 }
1567 
1568 
1570 {
1571  checkCellSubset();
1572 
1573  return fvMeshSubsetPtr_();
1574 }
1575 
1576 
1578 {
1579  checkCellSubset();
1580 
1581  return fvMeshSubsetPtr_();
1582 }
1583 
1584 
1586 {
1587  checkCellSubset();
1588 
1589  return pointMap_;
1590 }
1591 
1592 
1594 {
1595  checkCellSubset();
1596 
1597  return faceMap_;
1598 }
1599 
1600 
1602 {
1603  if (!faceFlipMapPtr_.valid())
1604  {
1605  const labelList& subToBaseFace = faceMap();
1606  const labelList& subToBaseCell = cellMap();
1607 
1608  faceFlipMapPtr_.reset(new labelList(subToBaseFace.size()));
1609  labelList& faceFlipMap = faceFlipMapPtr_();
1610 
1611  // Only exposed internal faces might be flipped (since we don't do
1612  // any cell renumbering, just compacting)
1613  label subInt = subMesh().nInternalFaces();
1614  const labelList& subOwn = subMesh().faceOwner();
1615  const labelList& own = baseMesh_.faceOwner();
1616 
1617  for (label subFacei = 0; subFacei < subInt; subFacei++)
1618  {
1619  faceFlipMap[subFacei] = subToBaseFace[subFacei] + 1;
1620  }
1621  for (label subFacei = subInt; subFacei < subOwn.size(); subFacei++)
1622  {
1623  const label facei = subToBaseFace[subFacei];
1624  if (subToBaseCell[subOwn[subFacei]] == own[facei])
1625  {
1626  faceFlipMap[subFacei] = facei + 1;
1627  }
1628  else
1629  {
1630  faceFlipMap[subFacei] = -facei - 1;
1631  }
1632  }
1633  }
1634 
1635  return faceFlipMapPtr_();
1636 }
1637 
1638 
1640 {
1641  checkCellSubset();
1642 
1643  return cellMap_;
1644 }
1645 
1646 
1648 {
1649  checkCellSubset();
1650 
1651  return patchMap_;
1652 }
1653 
1654 
1655 // ************************************************************************* //
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:202
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:113
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:493
const labelList & cellMap() const
Return cell map.
fvMeshSubset(const fvMesh &)
Construct given a mesh to subset.
Definition: fvMeshSubset.C:472
bool hasSubMesh() const
Have subMesh?
const labelList & faceFlipMap() const
Return face map with sign to encode flipped faces.
virtual ~fvMeshSubset()
Destructor.
Definition: fvMeshSubset.C:486
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:896
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:101
Constraint patch to hold internal faces exposed by sub-setting.
Foam::polyBoundaryMesh.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
label findPatchID(const word &patchName) const
Find patch index given a name.
wordList names() const
Return a list of patch names.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H: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:179
labelList getExposedFaces(const labelList &cellsToRemove) const
Get labels of exposed faces.
Definition: removeCells.C:73
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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:251
errorManip< error > abort(error &err)
Definition: errorManip.H:131
MeshZones< cellZone, polyMesh > meshCellZones
A MeshZones with the type cellZone.
messageStream Info
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
MeshZones< pointZone, polyMesh > meshPointZones
A MeshZones with the type pointZone.
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:211
error FatalError
static const label labelMax
Definition: label.H:62
MeshZones< faceZone, polyMesh > meshFaceZones
A MeshZones with the type faceZone.
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
wordList patchNames(nPatches)