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-2018 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 "polyRemoveFace.H"
30 #include "polyAddFace.H"
31 #include "polyModifyFace.H"
32 #include "polyRemovePoint.H"
33 #include "polyAddPoint.H"
34 #include "syncTools.H"
35 #include "meshTools.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(combineFaces, 0);
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 // Test face for (almost) convexeness. Allows certain convexity before
48 // complaining.
49 // For every two consecutive edges calculate the normal. If it differs too
50 // much from the average face normal complain.
51 bool Foam::combineFaces::convexFace
52 (
53  const scalar minConcaveCos,
54  const pointField& points,
55  const face& f
56 )
57 {
58  // Get outwards pointing normal of f.
59  const vector n = f.normal(points);
60 
61  // Get edge from f[0] to f[size-1];
62  vector ePrev(points[f.first()] - points[f.last()]);
63  scalar magEPrev = mag(ePrev);
64  ePrev /= magEPrev + vSmall;
65 
66  forAll(f, fp0)
67  {
68  // Get vertex after fp
69  label fp1 = f.fcIndex(fp0);
70 
71  // Normalized vector between two consecutive points
72  vector e10(points[f[fp1]] - points[f[fp0]]);
73  scalar magE10 = mag(e10);
74  e10 /= magE10 + vSmall;
75 
76  if (magEPrev > small && magE10 > small)
77  {
78  vector edgeNormal = ePrev ^ e10;
79 
80  if ((edgeNormal & n) < 0)
81  {
82  // Concave. Check angle.
83  if ((ePrev & e10) < minConcaveCos)
84  {
85  return false;
86  }
87  }
88  }
89 
90  ePrev = e10;
91  magEPrev = magE10;
92  }
93 
94  // Not a single internal angle is concave so face is convex.
95  return true;
96 }
97 
98 
99 // Determines if set of faces is valid to collapse into single face.
100 bool Foam::combineFaces::validFace
101 (
102  const scalar minConcaveCos,
103  const indirectPrimitivePatch& bigFace
104 )
105 {
106  // Get outside vertices (in local vertex numbering)
107  const labelListList& edgeLoops = bigFace.edgeLoops();
108 
109  if (edgeLoops.size() > 1)
110  {
111  return false;
112  }
113 
114  bool isNonManifold = bigFace.checkPointManifold(false, nullptr);
115  if (isNonManifold)
116  {
117  return false;
118  }
119 
120  // Check for convexness
121  face f(getOutsideFace(bigFace));
122 
123  return convexFace(minConcaveCos, bigFace.points(), f);
124 }
125 
126 
127 void Foam::combineFaces::regioniseFaces
128 (
129  const scalar minCos,
130  const label celli,
131  const labelList& cEdges,
132  Map<label>& faceRegion
133 ) const
134 {
135  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
136 
137  forAll(cEdges, i)
138  {
139  label edgeI = cEdges[i];
140 
141  label f0, f1;
142  meshTools::getEdgeFaces(mesh_, celli, edgeI, f0, f1);
143 
144  label p0 = patches.whichPatch(f0);
145  label p1 = patches.whichPatch(f1);
146 
147  // Face can be merged if
148  // - same non-coupled patch
149  // - small angle
150  if (p0 != -1 && p0 == p1 && !patches[p0].coupled())
151  {
152  vector f0Normal = mesh_.faceAreas()[f0];
153  f0Normal /= mag(f0Normal);
154  vector f1Normal = mesh_.faceAreas()[f1];
155  f1Normal /= mag(f1Normal);
156 
157  if ((f0Normal&f1Normal) > minCos)
158  {
159  Map<label>::const_iterator f0Fnd = faceRegion.find(f0);
160 
161  label region0 = -1;
162  if (f0Fnd != faceRegion.end())
163  {
164  region0 = f0Fnd();
165  }
166 
167  Map<label>::const_iterator f1Fnd = faceRegion.find(f1);
168 
169  label region1 = -1;
170  if (f1Fnd != faceRegion.end())
171  {
172  region1 = f1Fnd();
173  }
174 
175  if (region0 == -1)
176  {
177  if (region1 == -1)
178  {
179  label useRegion = faceRegion.size();
180  faceRegion.insert(f0, useRegion);
181  faceRegion.insert(f1, useRegion);
182  }
183  else
184  {
185  faceRegion.insert(f0, region1);
186  }
187  }
188  else
189  {
190  if (region1 == -1)
191  {
192  faceRegion.insert(f1, region0);
193  }
194  else if (region0 != region1)
195  {
196  // Merge the two regions
197  label useRegion = min(region0, region1);
198  label freeRegion = max(region0, region1);
199 
200  forAllIter(Map<label>, faceRegion, iter)
201  {
202  if (iter() == freeRegion)
203  {
204  iter() = useRegion;
205  }
206  }
207  }
208  }
209  }
210  }
211  }
212 }
213 
214 
215 bool Foam::combineFaces::faceNeighboursValid
216 (
217  const label celli,
218  const Map<label>& faceRegion
219 ) const
220 {
221  if (faceRegion.size() <= 1)
222  {
223  return true;
224  }
225 
226  const cell& cFaces = mesh_.cells()[celli];
227 
228  DynamicList<label> storage;
229 
230  // Test for face collapsing to edge since too many neighbours merged.
231  forAll(cFaces, cFacei)
232  {
233  label facei = cFaces[cFacei];
234 
235  if (!faceRegion.found(facei))
236  {
237  const labelList& fEdges = mesh_.faceEdges(facei, storage);
238 
239  // Count number of remaining faces neighbouring facei. This has
240  // to be 3 or more.
241 
242  // Unregioned neighbouring faces
243  DynamicList<label> neighbourFaces(cFaces.size());
244  // Regioned neighbouring faces
245  labelHashSet neighbourRegions;
246 
247  forAll(fEdges, i)
248  {
249  label edgeI = fEdges[i];
250  label nbrI = meshTools::otherFace(mesh_, celli, facei, edgeI);
251 
252  Map<label>::const_iterator iter = faceRegion.find(nbrI);
253 
254  if (iter == faceRegion.end())
255  {
256  if (findIndex(neighbourFaces, nbrI) == -1)
257  {
258  neighbourFaces.append(nbrI);
259  }
260  }
261  else
262  {
263  neighbourRegions.insert(iter());
264  }
265  }
266 
267  if ((neighbourFaces.size()+neighbourRegions.size()) < 3)
268  {
269  return false;
270  }
271  }
272  }
273  return true;
274 }
275 
276 
277 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
278 
279 // Construct from mesh
280 Foam::combineFaces::combineFaces
281 (
282  const polyMesh& mesh,
283  const bool undoable
284 )
285 :
286  mesh_(mesh),
287  undoable_(undoable),
288  masterFace_(0),
289  faceSetsVertices_(0),
290  savedPointLabels_(0),
291  savedPoints_(0)
292 {}
293 
294 
295 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
296 
298 (
299  const scalar featureCos,
300  const scalar minConcaveCos,
301  const labelHashSet& boundaryCells
302 ) const
303 {
304  // Lists of faces that can be merged.
305  DynamicList<labelList> allFaceSets(boundaryCells.size() / 10);
306  // Storage for on-the-fly cell-edge addressing.
307  DynamicList<label> storage;
308 
309  // On all cells regionise the faces
310  forAllConstIter(labelHashSet, boundaryCells, iter)
311  {
312  label celli = iter.key();
313 
314  const cell& cFaces = mesh_.cells()[celli];
315 
316  const labelList& cEdges = mesh_.cellEdges(celli, storage);
317 
318  // Region per face
319  Map<label> faceRegion(cFaces.size());
320  regioniseFaces(featureCos, celli, cEdges, faceRegion);
321 
322  // Now we have in faceRegion for every face the region with planar
323  // face sharing the same region. We now check whether the resulting
324  // sets cause a face
325  // - to become a set of edges since too many faces are merged.
326  // - to become convex
327 
328  if (faceNeighboursValid(celli, faceRegion))
329  {
330  // Create region-to-faces addressing
331  Map<labelList> regionToFaces(faceRegion.size());
332 
333  forAllConstIter(Map<label>, faceRegion, iter)
334  {
335  label facei = iter.key();
336  label region = iter();
337 
338  Map<labelList>::iterator regionFnd = regionToFaces.find(region);
339 
340  if (regionFnd != regionToFaces.end())
341  {
342  labelList& setFaces = regionFnd();
343  label sz = setFaces.size();
344  setFaces.setSize(sz+1);
345  setFaces[sz] = facei;
346  }
347  else
348  {
349  regionToFaces.insert(region, labelList(1, facei));
350  }
351  }
352 
353  // For every set check if it forms a valid convex face
354 
355  forAllConstIter(Map<labelList>, regionToFaces, iter)
356  {
357  // Make face out of setFaces
358  indirectPrimitivePatch bigFace
359  (
361  (
362  mesh_.faces(),
363  iter()
364  ),
365  mesh_.points()
366  );
367 
368  // Only store if -only one outside loop -which forms convex face
369  if (validFace(minConcaveCos, bigFace))
370  {
371  allFaceSets.append(iter());
372  }
373  }
374  }
375  }
376 
377  return allFaceSets.shrink();
378 }
379 
380 
382 (
383  const scalar featureCos,
384  const scalar minConcaveCos
385 ) const
386 {
387  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
388 
389  // Pick up all cells on boundary
390  labelHashSet boundaryCells(mesh_.nFaces()-mesh_.nInternalFaces());
391 
392  forAll(patches, patchi)
393  {
394  const polyPatch& patch = patches[patchi];
395 
396  if (!patch.coupled())
397  {
398  forAll(patch, i)
399  {
400  boundaryCells.insert(mesh_.faceOwner()[patch.start()+i]);
401  }
402  }
403  }
404 
405  return getMergeSets(featureCos, minConcaveCos, boundaryCells);
406 }
407 
408 
409 // Gets outside edgeloop as a face
410 // - in same order as faces
411 // - in mesh vertex labels
413 (
414  const indirectPrimitivePatch& fp
415 )
416 {
417  if (fp.edgeLoops().size() != 1)
418  {
420  << "Multiple outside loops:" << fp.edgeLoops()
421  << abort(FatalError);
422  }
423 
424  // Get first boundary edge. Since guaranteed one edgeLoop when in here this
425  // edge must be on it.
426  label bEdgeI = fp.nInternalEdges();
427 
428  const edge& e = fp.edges()[bEdgeI];
429 
430  const labelList& eFaces = fp.edgeFaces()[bEdgeI];
431 
432  if (eFaces.size() != 1)
433  {
435  << "boundary edge:" << bEdgeI
436  << " points:" << fp.meshPoints()[e[0]]
437  << ' ' << fp.meshPoints()[e[1]]
438  << " on indirectPrimitivePatch has " << eFaces.size()
439  << " faces using it" << abort(FatalError);
440  }
441 
442 
443  // Outside loop
444  const labelList& outsideLoop = fp.edgeLoops()[0];
445 
446 
447  // Get order of edge e in outside loop
448  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
449 
450  // edgeLoopConsistent : true = edge in same order as outsideloop
451  // false = opposite order
452  bool edgeLoopConsistent = false;
453 
454  {
455  label index0 = findIndex(outsideLoop, e[0]);
456  label index1 = findIndex(outsideLoop, e[1]);
457 
458  if (index0 == -1 || index1 == -1)
459  {
461  << "Cannot find boundary edge:" << e
462  << " points:" << fp.meshPoints()[e[0]]
463  << ' ' << fp.meshPoints()[e[1]]
464  << " in edgeLoop:" << outsideLoop << abort(FatalError);
465  }
466  else if (index1 == outsideLoop.fcIndex(index0))
467  {
468  edgeLoopConsistent = true;
469  }
470  else if (index0 == outsideLoop.fcIndex(index1))
471  {
472  edgeLoopConsistent = false;
473  }
474  else
475  {
477  << "Cannot find boundary edge:" << e
478  << " points:" << fp.meshPoints()[e[0]]
479  << ' ' << fp.meshPoints()[e[1]]
480  << " on consecutive points in edgeLoop:"
481  << outsideLoop << abort(FatalError);
482  }
483  }
484 
485 
486  // Get face in local vertices
487  const face& localF = fp.localFaces()[eFaces[0]];
488 
489  // Get order of edge in localF
490  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
491 
492  // faceEdgeConsistent : true = edge in same order as localF
493  // false = opposite order
494  bool faceEdgeConsistent = false;
495 
496  {
497  // Find edge in face.
498  label index = findIndex(fp.faceEdges()[eFaces[0]], bEdgeI);
499 
500  if (index == -1)
501  {
503  << "Cannot find boundary edge:" << e
504  << " points:" << fp.meshPoints()[e[0]]
505  << ' ' << fp.meshPoints()[e[1]]
506  << " in face:" << eFaces[0]
507  << " edges:" << fp.faceEdges()[eFaces[0]]
508  << abort(FatalError);
509  }
510  else if (localF[index] == e[0] && localF.nextLabel(index) == e[1])
511  {
512  faceEdgeConsistent = true;
513  }
514  else if (localF[index] == e[1] && localF.nextLabel(index) == e[0])
515  {
516  faceEdgeConsistent = false;
517  }
518  else
519  {
521  << "Cannot find boundary edge:" << e
522  << " points:" << fp.meshPoints()[e[0]]
523  << ' ' << fp.meshPoints()[e[1]]
524  << " in face:" << eFaces[0] << " verts:" << localF
525  << abort(FatalError);
526  }
527  }
528 
529  // Get face in mesh points.
530  face meshFace(renumber(fp.meshPoints(), outsideLoop));
531 
532  if (faceEdgeConsistent != edgeLoopConsistent)
533  {
534  reverse(meshFace);
535  }
536  return meshFace;
537 }
538 
539 
541 (
542  const labelListList& faceSets,
543  polyTopoChange& meshMod
544 )
545 {
546  if (undoable_)
547  {
548  masterFace_.setSize(faceSets.size());
549  faceSetsVertices_.setSize(faceSets.size());
550  forAll(faceSets, setI)
551  {
552  const labelList& setFaces = faceSets[setI];
553 
554  masterFace_[setI] = setFaces[0];
555  faceSetsVertices_[setI] = UIndirectList<face>
556  (
557  mesh_.faces(),
558  setFaces
559  );
560  }
561  }
562 
563  // Running count of number of faces using a point
564  labelList nPointFaces(mesh_.nPoints(), 0);
565 
566  const labelListList& pointFaces = mesh_.pointFaces();
567 
568  forAll(pointFaces, pointi)
569  {
570  nPointFaces[pointi] = pointFaces[pointi].size();
571  }
572 
573  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
574 
575  forAll(faceSets, setI)
576  {
577  const labelList& setFaces = faceSets[setI];
578 
579  // Check
580  if (debug)
581  {
582  forAll(setFaces, i)
583  {
584  label patchi = patches.whichPatch(setFaces[i]);
585 
586  if (patchi == -1 || patches[patchi].coupled())
587  {
589  << "Can only merge non-coupled boundary faces"
590  << " but found internal or coupled face:"
591  << setFaces[i] << " in set " << setI
592  << abort(FatalError);
593  }
594  }
595  }
596 
597  // Make face out of setFaces
598  indirectPrimitivePatch bigFace
599  (
601  (
602  mesh_.faces(),
603  setFaces
604  ),
605  mesh_.points()
606  );
607 
608  // Get outside vertices (in local vertex numbering)
609  const labelListList& edgeLoops = bigFace.edgeLoops();
610 
611  if (edgeLoops.size() != 1)
612  {
614  << "Faces to-be-merged " << setFaces
615  << " do not form a single big face." << nl
616  << abort(FatalError);
617  }
618 
619 
620  // Delete all faces except master, whose face gets modified.
621 
622  // Modify master face
623  // ~~~~~~~~~~~~~~~~~~
624 
625  label masterFacei = setFaces[0];
626 
627  // Get outside face in mesh vertex labels
628  face outsideFace(getOutsideFace(bigFace));
629 
630  label zoneID = mesh_.faceZones().whichZone(masterFacei);
631 
632  bool zoneFlip = false;
633 
634  if (zoneID >= 0)
635  {
636  const faceZone& fZone = mesh_.faceZones()[zoneID];
637 
638  zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
639  }
640 
641  label patchi = mesh_.boundaryMesh().whichPatch(masterFacei);
642 
643  meshMod.setAction
644  (
646  (
647  outsideFace, // modified face
648  masterFacei, // label of face being modified
649  mesh_.faceOwner()[masterFacei], // owner
650  -1, // neighbour
651  false, // face flip
652  patchi, // patch for face
653  false, // remove from zone
654  zoneID, // zone for face
655  zoneFlip // face flip in zone
656  )
657  );
658 
659 
660  // Delete all non-master faces
661  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
662 
663  for (label i = 1; i < setFaces.size(); i++)
664  {
665  meshMod.setAction(polyRemoveFace(setFaces[i]));
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.setAction(polyRemovePoint(pointi));
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.setAction(polyRemovePoint(pointi));
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.setAction
891  (
893  (
894  savedPoints_[localI], // point
895  -1, // master point
896  -1, // zone for point
897  true // supports a cell
898  )
899  );
900  restoredPoints.insert
901  (
902  addedPoints[localI], // current point label
903  savedPointLabels_[localI] // point label when it
904  // was stored
905  );
906  }
907  f[fp] = addedPoints[localI];
908  }
909  }
910  }
911 
912 
913  // Restore
914  // ~~~~~~~
915 
916  label own = mesh_.faceOwner()[masterFacei];
917  label zoneID = mesh_.faceZones().whichZone(masterFacei);
918  bool zoneFlip = false;
919  if (zoneID >= 0)
920  {
921  const faceZone& fZone = mesh_.faceZones()[zoneID];
922  zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
923  }
924  label patchi = mesh_.boundaryMesh().whichPatch(masterFacei);
925 
926  if (mesh_.boundaryMesh()[patchi].coupled())
927  {
929  << "Master face " << masterFacei << " is on coupled patch "
930  << mesh_.boundaryMesh()[patchi].name()
931  << abort(FatalError);
932  }
933 
934  // Pout<< "Restoring new master face " << masterFacei
935  // << " to vertices " << faces[0] << endl;
936 
937  // Modify the master face.
938  meshMod.setAction
939  (
941  (
942  faces[0], // original face
943  masterFacei, // label of face
944  own, // owner
945  -1, // neighbour
946  false, // face flip
947  patchi, // patch for face
948  false, // remove from zone
949  zoneID, // zone for face
950  zoneFlip // face flip in zone
951  )
952  );
953  restoredFaces.insert(masterFacei, masterFacei);
954 
955  // Add the previously removed faces
956  for (label i = 1; i < faces.size(); i++)
957  {
958  // Pout<< "Restoring removed face with vertices " << faces[i]
959  // << endl;
960 
961  label facei = meshMod.setAction
962  (
964  (
965  faces[i], // vertices
966  own, // owner,
967  -1, // neighbour,
968  -1, // masterPointID,
969  -1, // masterEdgeID,
970  masterFacei, // masterFaceID,
971  false, // flipFaceFlux,
972  patchi, // patchID,
973  zoneID, // zoneID,
974  zoneFlip // zoneFlip
975  )
976  );
977  restoredFaces.insert(facei, masterFacei);
978  }
979 
980  // Clear out restored set
981  faceSetsVertices_[setI].clear();
982  masterFace_[setI] = -1;
983  }
984 }
985 
986 
987 // ************************************************************************* //
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
Class containing data for face removal.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:463
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
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 inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Class describing modification of a face.
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:110
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:253
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:117
label nInternalEdges() const
Number of internal edges.
HashTable< label, label, Hash< label > >::const_iterator const_iterator
Definition: Map.H:59
Class containing data for point addition.
Definition: polyAddPoint.H:47
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:306
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
PrimitivePatch< face, IndirectList, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
scalar f1
Definition: createFields.H:28
const labelListList & edgeLoops() const
Return list of closed loops of boundary vertices.
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:51
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
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:142
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
An STL-conforming iterator.
Definition: HashTable.H:426
static face getOutsideFace(const indirectPrimitivePatch &)
Gets outside of patch as a face (in mesh point labels)
Definition: combineFaces.C:413
A list of faces which address into the list of points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:310
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
label size() const
Return the number of elements in the list.
const labelListList & edgeFaces() const
Return edge-face addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::polyBoundaryMesh.
labelListList getMergeSets(const scalar featureCos, const scalar minConcaveCos, const labelHashSet &boundaryCells) const
Extract lists of all (non-coupled) boundary faces on selected.
Definition: combineFaces.C:298
void setRefinement(const labelListList &, polyTopoChange &)
Play commands into polyTopoChange to combine faces. Gets.
Definition: combineFaces.C:541
void reverse(UList< T > &, const label n)
Definition: UListI.H:322
static const char nl
Definition: Ostream.H:265
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
A List with indirect addressing.
Definition: fvMatrix.H:106
Direct mesh changes based on v1.3 polyTopoChange syntax.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
const labelListList & faceEdges() const
Return face-edge addressing.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: combineFaces.C:766
A List with indirect addressing.
Definition: IndirectList.H:102
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
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
Class containing data for point removal.
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
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:495
Namespace for OpenFOAM.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.