combineFaces.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 \*---------------------------------------------------------------------------*/
25 
26 #include "combineFaces.H"
27 #include "polyMesh.H"
28 #include "polyTopoChange.H"
29 #include "syncTools.H"
30 #include "meshTools.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 // Test face for (almost) convexeness. Allows certain convexity before
43 // complaining.
44 // For every two consecutive edges calculate the normal. If it differs too
45 // much from the average face normal complain.
46 bool Foam::combineFaces::convexFace
47 (
48  const scalar minConcaveCos,
49  const pointField& points,
50  const face& f
51 )
52 {
53  // Get outwards pointing normal of f.
54  const vector n = f.normal(points);
55 
56  // Get edge from f[0] to f[size-1];
57  vector ePrev(points[f.first()] - points[f.last()]);
58  scalar magEPrev = mag(ePrev);
59  ePrev /= magEPrev + vSmall;
60 
61  forAll(f, fp0)
62  {
63  // Get vertex after fp
64  label fp1 = f.fcIndex(fp0);
65 
66  // Normalised vector between two consecutive points
67  vector e10(points[f[fp1]] - points[f[fp0]]);
68  scalar magE10 = mag(e10);
69  e10 /= magE10 + vSmall;
70 
71  if (magEPrev > small && magE10 > small)
72  {
73  vector edgeNormal = ePrev ^ e10;
74 
75  if ((edgeNormal & n) < 0)
76  {
77  // Concave. Check angle.
78  if ((ePrev & e10) < minConcaveCos)
79  {
80  return false;
81  }
82  }
83  }
84 
85  ePrev = e10;
86  magEPrev = magE10;
87  }
88 
89  // Not a single internal angle is concave so face is convex.
90  return true;
91 }
92 
93 
94 // Determines if set of faces is valid to collapse into single face.
95 bool Foam::combineFaces::validFace
96 (
97  const scalar minConcaveCos,
98  const indirectPrimitivePatch& bigFace
99 )
100 {
101  // Get outside vertices (in local vertex numbering)
102  const labelListList& edgeLoops = bigFace.edgeLoops();
103 
104  if (edgeLoops.size() > 1)
105  {
106  return false;
107  }
108 
109  bool isNonManifold = bigFace.checkPointManifold(false, nullptr);
110  if (isNonManifold)
111  {
112  return false;
113  }
114 
115  // Check for convexness
116  face f(getOutsideFace(bigFace));
117 
118  return convexFace(minConcaveCos, bigFace.points(), f);
119 }
120 
121 
122 void Foam::combineFaces::regioniseFaces
123 (
124  const scalar minCos,
125  const labelHashSet& patchIDs,
126  const label celli,
127  const labelList& cEdges,
128  Map<label>& faceRegion
129 ) const
130 {
131  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
132 
133  forAll(cEdges, i)
134  {
135  label edgeI = cEdges[i];
136 
137  label f0, f1;
138  meshTools::getEdgeFaces(mesh_, celli, edgeI, f0, f1);
139 
140  label p0 = patches.whichPatch(f0);
141  label p1 = patches.whichPatch(f1);
142 
143  // Face can be merged if
144  // - same non-coupled patch in list
145  // - small angle
146  if
147  (
148  p0 != -1
149  && p0 == p1
150  && !patches[p0].coupled()
151  && patchIDs.found(p0)
152  )
153  {
154  vector f0Normal = mesh_.faceAreas()[f0];
155  f0Normal /= mag(f0Normal);
156  vector f1Normal = mesh_.faceAreas()[f1];
157  f1Normal /= mag(f1Normal);
158 
159  if ((f0Normal&f1Normal) > minCos)
160  {
161  Map<label>::const_iterator f0Fnd = faceRegion.find(f0);
162 
163  label region0 = -1;
164  if (f0Fnd != faceRegion.end())
165  {
166  region0 = f0Fnd();
167  }
168 
169  Map<label>::const_iterator f1Fnd = faceRegion.find(f1);
170 
171  label region1 = -1;
172  if (f1Fnd != faceRegion.end())
173  {
174  region1 = f1Fnd();
175  }
176 
177  if (region0 == -1)
178  {
179  if (region1 == -1)
180  {
181  label useRegion = faceRegion.size();
182  faceRegion.insert(f0, useRegion);
183  faceRegion.insert(f1, useRegion);
184  }
185  else
186  {
187  faceRegion.insert(f0, region1);
188  }
189  }
190  else
191  {
192  if (region1 == -1)
193  {
194  faceRegion.insert(f1, region0);
195  }
196  else if (region0 != region1)
197  {
198  // Merge the two regions
199  label useRegion = min(region0, region1);
200  label freeRegion = max(region0, region1);
201 
202  forAllIter(Map<label>, faceRegion, iter)
203  {
204  if (iter() == freeRegion)
205  {
206  iter() = useRegion;
207  }
208  }
209  }
210  }
211  }
212  }
213  }
214 }
215 
216 
217 bool Foam::combineFaces::faceNeighboursValid
218 (
219  const label celli,
220  const Map<label>& faceRegion
221 ) const
222 {
223  if (faceRegion.size() <= 1)
224  {
225  return true;
226  }
227 
228  const cell& cFaces = mesh_.cells()[celli];
229 
230  DynamicList<label> storage;
231 
232  // Test for face collapsing to edge since too many neighbours merged.
233  forAll(cFaces, cFacei)
234  {
235  label facei = cFaces[cFacei];
236 
237  if (!faceRegion.found(facei))
238  {
239  const labelList& fEdges = mesh_.faceEdges(facei, storage);
240 
241  // Count number of remaining faces neighbouring facei. This has
242  // to be 3 or more.
243 
244  // Unregioned neighbouring faces
245  DynamicList<label> neighbourFaces(cFaces.size());
246  // Regioned neighbouring faces
247  labelHashSet neighbourRegions;
248 
249  forAll(fEdges, i)
250  {
251  label edgeI = fEdges[i];
252  label nbrI = meshTools::otherFace(mesh_, celli, facei, edgeI);
253 
254  Map<label>::const_iterator iter = faceRegion.find(nbrI);
255 
256  if (iter == faceRegion.end())
257  {
258  if (findIndex(neighbourFaces, nbrI) == -1)
259  {
260  neighbourFaces.append(nbrI);
261  }
262  }
263  else
264  {
265  neighbourRegions.insert(iter());
266  }
267  }
268 
269  if ((neighbourFaces.size()+neighbourRegions.size()) < 3)
270  {
271  return false;
272  }
273  }
274  }
275  return true;
276 }
277 
278 
279 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
280 
281 // Construct from mesh
283 (
284  const polyMesh& mesh,
285  const bool undoable
286 )
287 :
288  mesh_(mesh),
289  undoable_(undoable),
290  masterFace_(0),
291  faceSetsVertices_(0),
292  savedPointLabels_(0),
293  savedPoints_(0)
294 {}
295 
296 
297 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
298 
300 (
301  const scalar featureCos,
302  const scalar minConcaveCos,
303  const labelHashSet& patchIDs,
304  const labelHashSet& boundaryCells
305 ) const
306 {
307  // Lists of faces that can be merged.
308  DynamicList<labelList> allFaceSets(boundaryCells.size() / 10);
309  // Storage for on-the-fly cell-edge addressing.
310  DynamicList<label> storage;
311 
312  // On all cells regionise the faces
313  forAllConstIter(labelHashSet, boundaryCells, iter)
314  {
315  label celli = iter.key();
316 
317  const cell& cFaces = mesh_.cells()[celli];
318 
319  const labelList& cEdges = mesh_.cellEdges(celli, storage);
320 
321  // Region per face
322  Map<label> faceRegion(cFaces.size());
323  regioniseFaces(featureCos, patchIDs, celli, cEdges, faceRegion);
324 
325  // Now we have in faceRegion for every face the region with planar
326  // face sharing the same region. We now check whether the resulting
327  // sets cause a face
328  // - to become a set of edges since too many faces are merged.
329  // - to become convex
330 
331  if (faceNeighboursValid(celli, faceRegion))
332  {
333  // Create region-to-faces addressing
334  Map<labelList> regionToFaces(faceRegion.size());
335 
336  forAllConstIter(Map<label>, faceRegion, iter)
337  {
338  label facei = iter.key();
339  label region = iter();
340 
341  Map<labelList>::iterator regionFnd = regionToFaces.find(region);
342 
343  if (regionFnd != regionToFaces.end())
344  {
345  labelList& setFaces = regionFnd();
346  label sz = setFaces.size();
347  setFaces.setSize(sz+1);
348  setFaces[sz] = facei;
349  }
350  else
351  {
352  regionToFaces.insert(region, labelList(1, facei));
353  }
354  }
355 
356  // For every set check if it forms a valid convex face
357 
358  forAllConstIter(Map<labelList>, regionToFaces, iter)
359  {
360  // Make face out of setFaces
361  indirectPrimitivePatch bigFace
362  (
364  (
365  mesh_.faces(),
366  iter()
367  ),
368  mesh_.points()
369  );
370 
371  // Only store if -only one outside loop -which forms convex face
372  if (validFace(minConcaveCos, bigFace))
373  {
374  allFaceSets.append(iter());
375  }
376  }
377  }
378  }
379 
380  return allFaceSets.shrink();
381 }
382 
383 
385 (
386  const scalar featureCos,
387  const scalar minConcaveCos,
388  const labelHashSet& patchIDs
389 ) const
390 {
391  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
392 
393  // Pick up all cells on boundary
394  labelHashSet boundaryCells(mesh_.nFaces()-mesh_.nInternalFaces());
395 
396  forAllConstIter(labelHashSet, patchIDs, iter)
397  {
398  const polyPatch& patch = patches[iter.key()];
399 
400  if (!patch.coupled())
401  {
402  forAll(patch, i)
403  {
404  boundaryCells.insert(mesh_.faceOwner()[patch.start()+i]);
405  }
406  }
407  }
408 
409  return getMergeSets(featureCos, minConcaveCos, patchIDs, boundaryCells);
410 }
411 
412 
413 
415 (
416  const scalar featureCos,
417  const scalar minConcaveCos
418 ) const
419 {
420  const labelHashSet patchIDs(identityMap(mesh_.boundaryMesh().size()));
421 
422  return getMergeSets(featureCos, minConcaveCos, patchIDs);
423 }
424 
425 
426 // Gets outside edgeloop as a face
427 // - in same order as faces
428 // - in mesh vertex labels
430 (
431  const indirectPrimitivePatch& fp
432 )
433 {
434  if (fp.edgeLoops().size() != 1)
435  {
437  << "Multiple outside loops:" << fp.edgeLoops()
438  << abort(FatalError);
439  }
440 
441  // Get first boundary edge. Since guaranteed one edgeLoop when in here this
442  // edge must be on it.
443  label bEdgeI = fp.nInternalEdges();
444 
445  const edge& e = fp.edges()[bEdgeI];
446 
447  const labelList& eFaces = fp.edgeFaces()[bEdgeI];
448 
449  if (eFaces.size() != 1)
450  {
452  << "boundary edge:" << bEdgeI
453  << " points:" << fp.meshPoints()[e[0]]
454  << ' ' << fp.meshPoints()[e[1]]
455  << " on indirectPrimitivePatch has " << eFaces.size()
456  << " faces using it" << abort(FatalError);
457  }
458 
459 
460  // Outside loop
461  const labelList& outsideLoop = fp.edgeLoops()[0];
462 
463 
464  // Get order of edge e in outside loop
465  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
466 
467  // edgeLoopConsistent : true = edge in same order as outsideloop
468  // false = opposite order
469  bool edgeLoopConsistent = false;
470 
471  {
472  label index0 = findIndex(outsideLoop, e[0]);
473  label index1 = findIndex(outsideLoop, e[1]);
474 
475  if (index0 == -1 || index1 == -1)
476  {
478  << "Cannot find boundary edge:" << e
479  << " points:" << fp.meshPoints()[e[0]]
480  << ' ' << fp.meshPoints()[e[1]]
481  << " in edgeLoop:" << outsideLoop << abort(FatalError);
482  }
483  else if (index1 == outsideLoop.fcIndex(index0))
484  {
485  edgeLoopConsistent = true;
486  }
487  else if (index0 == outsideLoop.fcIndex(index1))
488  {
489  edgeLoopConsistent = false;
490  }
491  else
492  {
494  << "Cannot find boundary edge:" << e
495  << " points:" << fp.meshPoints()[e[0]]
496  << ' ' << fp.meshPoints()[e[1]]
497  << " on consecutive points in edgeLoop:"
498  << outsideLoop << abort(FatalError);
499  }
500  }
501 
502 
503  // Get face in local vertices
504  const face& localF = fp.localFaces()[eFaces[0]];
505 
506  // Get order of edge in localF
507  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
508 
509  // faceEdgeConsistent : true = edge in same order as localF
510  // false = opposite order
511  bool faceEdgeConsistent = false;
512 
513  {
514  // Find edge in face.
515  label index = findIndex(fp.faceEdges()[eFaces[0]], bEdgeI);
516 
517  if (index == -1)
518  {
520  << "Cannot find boundary edge:" << e
521  << " points:" << fp.meshPoints()[e[0]]
522  << ' ' << fp.meshPoints()[e[1]]
523  << " in face:" << eFaces[0]
524  << " edges:" << fp.faceEdges()[eFaces[0]]
525  << abort(FatalError);
526  }
527  else if (localF[index] == e[0] && localF.nextLabel(index) == e[1])
528  {
529  faceEdgeConsistent = true;
530  }
531  else if (localF[index] == e[1] && localF.nextLabel(index) == e[0])
532  {
533  faceEdgeConsistent = false;
534  }
535  else
536  {
538  << "Cannot find boundary edge:" << e
539  << " points:" << fp.meshPoints()[e[0]]
540  << ' ' << fp.meshPoints()[e[1]]
541  << " in face:" << eFaces[0] << " verts:" << localF
542  << abort(FatalError);
543  }
544  }
545 
546  // Get face in mesh points.
547  face meshFace(renumber(fp.meshPoints(), outsideLoop));
548 
549  if (faceEdgeConsistent != edgeLoopConsistent)
550  {
551  reverse(meshFace);
552  }
553  return meshFace;
554 }
555 
556 
558 (
559  const labelListList& faceSets,
560  polyTopoChange& meshMod
561 )
562 {
563  if (undoable_)
564  {
565  masterFace_.setSize(faceSets.size());
566  faceSetsVertices_.setSize(faceSets.size());
567  forAll(faceSets, setI)
568  {
569  const labelList& setFaces = faceSets[setI];
570 
571  masterFace_[setI] = setFaces[0];
572  faceSetsVertices_[setI] = UIndirectList<face>
573  (
574  mesh_.faces(),
575  setFaces
576  );
577  }
578  }
579 
580  // Running count of number of faces using a point
581  labelList nPointFaces(mesh_.nPoints(), 0);
582 
583  const labelListList& pointFaces = mesh_.pointFaces();
584 
585  forAll(pointFaces, pointi)
586  {
587  nPointFaces[pointi] = pointFaces[pointi].size();
588  }
589 
590  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
591 
592  forAll(faceSets, setI)
593  {
594  const labelList& setFaces = faceSets[setI];
595 
596  // Check
597  if (debug)
598  {
599  forAll(setFaces, i)
600  {
601  label patchi = patches.whichPatch(setFaces[i]);
602 
603  if (patchi == -1 || patches[patchi].coupled())
604  {
606  << "Can only merge non-coupled boundary faces"
607  << " but found internal or coupled face:"
608  << setFaces[i] << " in set " << setI
609  << abort(FatalError);
610  }
611  }
612  }
613 
614  // Make face out of setFaces
615  indirectPrimitivePatch bigFace
616  (
618  (
619  mesh_.faces(),
620  setFaces
621  ),
622  mesh_.points()
623  );
624 
625  // Get outside vertices (in local vertex numbering)
626  const labelListList& edgeLoops = bigFace.edgeLoops();
627 
628  if (edgeLoops.size() != 1)
629  {
631  << "Faces to-be-merged " << setFaces
632  << " do not form a single big face." << nl
633  << abort(FatalError);
634  }
635 
636 
637  // Delete all faces except master, whose face gets modified.
638 
639  // Modify master face
640  // ~~~~~~~~~~~~~~~~~~
641 
642  const label masterFacei = setFaces[0];
643 
644  // Get outside face in mesh vertex labels
645  const face outsideFace(getOutsideFace(bigFace));
646 
647  const label patchi = mesh_.boundaryMesh().whichPatch(masterFacei);
648 
649  meshMod.modifyFace
650  (
651  outsideFace, // modified face
652  masterFacei, // label of face being modified
653  mesh_.faceOwner()[masterFacei], // owner
654  -1, // neighbour
655  false, // face flip
656  patchi // patch for face
657  );
658 
659 
660  // Delete all non-master faces
661  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
662 
663  for (label i = 1; i < setFaces.size(); i++)
664  {
665  meshMod.removeFace(setFaces[i], -1);
666  }
667 
668 
669  // Mark unused points for deletion
670  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
671 
672  // Remove count of all faces combined
673  forAll(setFaces, i)
674  {
675  const face& f = mesh_.faces()[setFaces[i]];
676 
677  forAll(f, fp)
678  {
679  nPointFaces[f[fp]]--;
680  }
681  }
682  // But keep count on new outside face
683  forAll(outsideFace, fp)
684  {
685  nPointFaces[outsideFace[fp]]++;
686  }
687  }
688 
689 
690  // If one or more processors use the point it needs to be kept.
692  (
693  mesh_,
694  nPointFaces,
695  plusEqOp<label>(),
696  label(0) // null value
697  );
698 
699  // Remove all unused points. Store position if undoable.
700  if (!undoable_)
701  {
702  forAll(nPointFaces, pointi)
703  {
704  if (nPointFaces[pointi] == 0)
705  {
706  meshMod.removePoint(pointi, -1);
707  }
708  }
709  }
710  else
711  {
712  // Count removed points
713  label n = 0;
714  forAll(nPointFaces, pointi)
715  {
716  if (nPointFaces[pointi] == 0)
717  {
718  n++;
719  }
720  }
721 
722  savedPointLabels_.setSize(n);
723  savedPoints_.setSize(n);
724  Map<label> meshToSaved(2*n);
725 
726  // Remove points and store position
727  n = 0;
728  forAll(nPointFaces, pointi)
729  {
730  if (nPointFaces[pointi] == 0)
731  {
732  meshMod.removePoint(pointi, -1);
733 
734  savedPointLabels_[n] = pointi;
735  savedPoints_[n] = mesh_.points()[pointi];
736 
737  meshToSaved.insert(pointi, n);
738  n++;
739  }
740  }
741 
742  // Update stored vertex labels. Negative indices index into local points
743  forAll(faceSetsVertices_, setI)
744  {
745  faceList& setFaces = faceSetsVertices_[setI];
746 
747  forAll(setFaces, i)
748  {
749  face& f = setFaces[i];
750 
751  forAll(f, fp)
752  {
753  label pointi = f[fp];
754 
755  if (nPointFaces[pointi] == 0)
756  {
757  f[fp] = -meshToSaved[pointi]-1;
758  }
759  }
760  }
761  }
762  }
763 }
764 
765 
767 {
768  if (undoable_)
769  {
770  // Master face just renumbering of point labels
771  inplaceRenumber(map.reverseFaceMap(), masterFace_);
772 
773  // Stored faces refer to backed-up vertices (not changed)
774  // and normal vertices (need to be renumbered)
775  forAll(faceSetsVertices_, setI)
776  {
777  faceList& faces = faceSetsVertices_[setI];
778 
779  forAll(faces, i)
780  {
781  // Note: faces[i] can have negative elements.
782  face& f = faces[i];
783 
784  forAll(f, fp)
785  {
786  label pointi = f[fp];
787 
788  if (pointi >= 0)
789  {
790  f[fp] = map.reversePointMap()[pointi];
791 
792  if (f[fp] < 0)
793  {
795  << "In set " << setI << " at position " << i
796  << " with master face "
797  << masterFace_[setI] << nl
798  << "the points of the slave face " << faces[i]
799  << " don't exist anymore."
800  << abort(FatalError);
801  }
802  }
803  }
804  }
805  }
806  }
807 }
808 
809 
810 // Note that faces on coupled patches are never combined (or at least
811 // getMergeSets never picks them up. Thus the points on them will never get
812 // deleted since still used by the coupled faces. So never a risk of undoing
813 // things (faces/points) on coupled patches.
815 (
816  const labelList& masterFaces,
817  polyTopoChange& meshMod,
818  Map<label>& restoredPoints,
819  Map<label>& restoredFaces,
820  Map<label>& restoredCells
821 )
822 {
823  if (!undoable_)
824  {
826  << "Can only call setUnrefinement if constructed with"
827  << " unrefinement capability." << exit(FatalError);
828  }
829 
830 
831  // Restored points
832  labelList addedPoints(savedPoints_.size(), -1);
833 
834  // Invert set-to-master-face
835  Map<label> masterToSet(masterFace_.size());
836 
837  forAll(masterFace_, setI)
838  {
839  if (masterFace_[setI] >= 0)
840  {
841  masterToSet.insert(masterFace_[setI], setI);
842  }
843  }
844 
845  forAll(masterFaces, i)
846  {
847  label masterFacei = masterFaces[i];
848 
849  Map<label>::const_iterator iter = masterToSet.find(masterFacei);
850 
851  if (iter == masterToSet.end())
852  {
854  << "Master face " << masterFacei
855  << " is not the master of one of the merge sets"
856  << " or has already been merged"
857  << abort(FatalError);
858  }
859 
860  label setI = iter();
861 
862 
863  // Update faces of the merge setI for reintroduced vertices
864  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
865 
866  faceList& faces = faceSetsVertices_[setI];
867 
868  if (faces.empty())
869  {
871  << "Set " << setI << " with master face " << masterFacei
872  << " has already been merged." << abort(FatalError);
873  }
874 
875  forAll(faces, i)
876  {
877  face& f = faces[i];
878 
879  forAll(f, fp)
880  {
881  label pointi = f[fp];
882 
883  if (pointi < 0)
884  {
885  label localI = -pointi-1;
886 
887  if (addedPoints[localI] == -1)
888  {
889  // First occurrence of saved point. Reintroduce point
890  addedPoints[localI] = meshMod.addPoint
891  (
892  savedPoints_[localI], // point
893  -1, // master point
894  true // supports a cell
895  );
896  restoredPoints.insert
897  (
898  addedPoints[localI], // current point label
899  savedPointLabels_[localI] // point label when it
900  // was stored
901  );
902  }
903  f[fp] = addedPoints[localI];
904  }
905  }
906  }
907 
908 
909  // Restore
910  // ~~~~~~~
911 
912  const label own = mesh_.faceOwner()[masterFacei];
913  const label patchi = mesh_.boundaryMesh().whichPatch(masterFacei);
914 
915  if (mesh_.boundaryMesh()[patchi].coupled())
916  {
918  << "Master face " << masterFacei << " is on coupled patch "
919  << mesh_.boundaryMesh()[patchi].name()
920  << abort(FatalError);
921  }
922 
923  // Pout<< "Restoring new master face " << masterFacei
924  // << " to vertices " << faces[0] << endl;
925 
926  // Modify the master face.
927  meshMod.modifyFace
928  (
929  faces[0], // original face
930  masterFacei, // label of face
931  own, // owner
932  -1, // neighbour
933  false, // face flip
934  patchi // patch for face
935  );
936  restoredFaces.insert(masterFacei, masterFacei);
937 
938  // Add the previously removed faces
939  for (label i = 1; i < faces.size(); i++)
940  {
941  // Pout<< "Restoring removed face with vertices " << faces[i]
942  // << endl;
943 
944  label facei = meshMod.addFace
945  (
946  faces[i], // vertices
947  own, // owner,
948  -1, // neighbour,
949  masterFacei, // masterFaceID,
950  false, // flipFaceFlux,
951  patchi // patchID,
952  );
953  restoredFaces.insert(facei, masterFacei);
954  }
955 
956  // Clear out restored set
957  faceSetsVertices_[setI].clear();
958  masterFace_[setI] = -1;
959  }
960 }
961 
962 
963 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
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
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:167
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:197
A List with indirect addressing.
Definition: IndirectList.H:105
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
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const labelListList & edgeLoops() const
Return list of closed loops of boundary vertices.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
const labelListList & edgeFaces() const
Return edge-face addressing.
const labelListList & faceEdges() const
Return face-edge addressing.
A List with indirect addressing.
Definition: UIndirectList.H:60
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
Combines boundary faces into single face. The faces get the patch of the first face ('the master')
Definition: combineFaces.H:56
labelListList getMergeSets(const scalar featureCos, const scalar minConcaveCos, const labelHashSet &patchIDs, const labelHashSet &boundaryCells) const
Extract lists of all (non-coupled) boundary faces on selected.
Definition: combineFaces.C:300
void topoChange(const polyTopoChangeMap &)
Force recalculation of locally stored data on topological change.
Definition: combineFaces.C:766
void setUnrefinement(const labelList &masterFaces, polyTopoChange &meshMod, Map< label > &restoredPoints, Map< label > &restoredFaces, Map< label > &restoredCells)
Play commands into polyTopoChange to reinsert original faces.
Definition: combineFaces.C:815
void setRefinement(const labelListList &, polyTopoChange &)
Play commands into polyTopoChange to combine faces. Gets.
Definition: combineFaces.C:558
combineFaces(const polyMesh &mesh, const bool undoable=false)
Construct from mesh.
Definition: combineFaces.C:283
static face getOutsideFace(const indirectPrimitivePatch &)
Gets outside of patch as a face (in mesh point labels)
Definition: combineFaces.C:430
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:100
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:290
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reversePointMap() const
Reverse point map.
const labelList & reverseFaceMap() const
Reverse face map.
Direct mesh changes based on v1.3 polyTopoChange syntax.
void removePoint(const label, const label)
Remove point / merge points.
void removeFace(const label, const label)
Remove face / merge faces.
label addPoint(const point &, const label masterPointID, const bool inCell)
Add point and return new point index.
void modifyFace(const face &f, const label facei, const label own, const label nei, const bool flipFaceFlux, const label patchID)
Modify vertices or cell of face.
label addFace(const face &f, const label own, const label nei, const label masterFaceID, const bool flipFaceFlux, const label patchID)
Add face to cells and return new face index.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const pointField & points
const fvPatchList & patches
void getEdgeFaces(const primitiveMesh &, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition: meshTools.C:458
label otherFace(const primitiveMesh &, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:536
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:106
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
static const char nl
Definition: Ostream.H:266
labelList f(nPoints)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112