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