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-2019 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 labelHashSet& patchIDs,
131  const label celli,
132  const labelList& cEdges,
133  Map<label>& faceRegion
134 ) const
135 {
136  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
137 
138  forAll(cEdges, i)
139  {
140  label edgeI = cEdges[i];
141 
142  label f0, f1;
143  meshTools::getEdgeFaces(mesh_, celli, edgeI, f0, f1);
144 
145  label p0 = patches.whichPatch(f0);
146  label p1 = patches.whichPatch(f1);
147 
148  // Face can be merged if
149  // - same non-coupled patch in list
150  // - small angle
151  if
152  (
153  p0 != -1
154  && p0 == p1
155  && !patches[p0].coupled()
156  && patchIDs.found(p0)
157  )
158  {
159  vector f0Normal = mesh_.faceAreas()[f0];
160  f0Normal /= mag(f0Normal);
161  vector f1Normal = mesh_.faceAreas()[f1];
162  f1Normal /= mag(f1Normal);
163 
164  if ((f0Normal&f1Normal) > minCos)
165  {
166  Map<label>::const_iterator f0Fnd = faceRegion.find(f0);
167 
168  label region0 = -1;
169  if (f0Fnd != faceRegion.end())
170  {
171  region0 = f0Fnd();
172  }
173 
174  Map<label>::const_iterator f1Fnd = faceRegion.find(f1);
175 
176  label region1 = -1;
177  if (f1Fnd != faceRegion.end())
178  {
179  region1 = f1Fnd();
180  }
181 
182  if (region0 == -1)
183  {
184  if (region1 == -1)
185  {
186  label useRegion = faceRegion.size();
187  faceRegion.insert(f0, useRegion);
188  faceRegion.insert(f1, useRegion);
189  }
190  else
191  {
192  faceRegion.insert(f0, region1);
193  }
194  }
195  else
196  {
197  if (region1 == -1)
198  {
199  faceRegion.insert(f1, region0);
200  }
201  else if (region0 != region1)
202  {
203  // Merge the two regions
204  label useRegion = min(region0, region1);
205  label freeRegion = max(region0, region1);
206 
207  forAllIter(Map<label>, faceRegion, iter)
208  {
209  if (iter() == freeRegion)
210  {
211  iter() = useRegion;
212  }
213  }
214  }
215  }
216  }
217  }
218  }
219 }
220 
221 
222 bool Foam::combineFaces::faceNeighboursValid
223 (
224  const label celli,
225  const Map<label>& faceRegion
226 ) const
227 {
228  if (faceRegion.size() <= 1)
229  {
230  return true;
231  }
232 
233  const cell& cFaces = mesh_.cells()[celli];
234 
235  DynamicList<label> storage;
236 
237  // Test for face collapsing to edge since too many neighbours merged.
238  forAll(cFaces, cFacei)
239  {
240  label facei = cFaces[cFacei];
241 
242  if (!faceRegion.found(facei))
243  {
244  const labelList& fEdges = mesh_.faceEdges(facei, storage);
245 
246  // Count number of remaining faces neighbouring facei. This has
247  // to be 3 or more.
248 
249  // Unregioned neighbouring faces
250  DynamicList<label> neighbourFaces(cFaces.size());
251  // Regioned neighbouring faces
252  labelHashSet neighbourRegions;
253 
254  forAll(fEdges, i)
255  {
256  label edgeI = fEdges[i];
257  label nbrI = meshTools::otherFace(mesh_, celli, facei, edgeI);
258 
259  Map<label>::const_iterator iter = faceRegion.find(nbrI);
260 
261  if (iter == faceRegion.end())
262  {
263  if (findIndex(neighbourFaces, nbrI) == -1)
264  {
265  neighbourFaces.append(nbrI);
266  }
267  }
268  else
269  {
270  neighbourRegions.insert(iter());
271  }
272  }
273 
274  if ((neighbourFaces.size()+neighbourRegions.size()) < 3)
275  {
276  return false;
277  }
278  }
279  }
280  return true;
281 }
282 
283 
284 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
285 
286 // Construct from mesh
288 (
289  const polyMesh& mesh,
290  const bool undoable
291 )
292 :
293  mesh_(mesh),
294  undoable_(undoable),
295  masterFace_(0),
296  faceSetsVertices_(0),
297  savedPointLabels_(0),
298  savedPoints_(0)
299 {}
300 
301 
302 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
303 
305 (
306  const scalar featureCos,
307  const scalar minConcaveCos,
308  const labelHashSet& patchIDs,
309  const labelHashSet& boundaryCells
310 ) const
311 {
312  // Lists of faces that can be merged.
313  DynamicList<labelList> allFaceSets(boundaryCells.size() / 10);
314  // Storage for on-the-fly cell-edge addressing.
315  DynamicList<label> storage;
316 
317  // On all cells regionise the faces
318  forAllConstIter(labelHashSet, boundaryCells, iter)
319  {
320  label celli = iter.key();
321 
322  const cell& cFaces = mesh_.cells()[celli];
323 
324  const labelList& cEdges = mesh_.cellEdges(celli, storage);
325 
326  // Region per face
327  Map<label> faceRegion(cFaces.size());
328  regioniseFaces(featureCos, patchIDs, celli, cEdges, faceRegion);
329 
330  // Now we have in faceRegion for every face the region with planar
331  // face sharing the same region. We now check whether the resulting
332  // sets cause a face
333  // - to become a set of edges since too many faces are merged.
334  // - to become convex
335 
336  if (faceNeighboursValid(celli, faceRegion))
337  {
338  // Create region-to-faces addressing
339  Map<labelList> regionToFaces(faceRegion.size());
340 
341  forAllConstIter(Map<label>, faceRegion, iter)
342  {
343  label facei = iter.key();
344  label region = iter();
345 
346  Map<labelList>::iterator regionFnd = regionToFaces.find(region);
347 
348  if (regionFnd != regionToFaces.end())
349  {
350  labelList& setFaces = regionFnd();
351  label sz = setFaces.size();
352  setFaces.setSize(sz+1);
353  setFaces[sz] = facei;
354  }
355  else
356  {
357  regionToFaces.insert(region, labelList(1, facei));
358  }
359  }
360 
361  // For every set check if it forms a valid convex face
362 
363  forAllConstIter(Map<labelList>, regionToFaces, iter)
364  {
365  // Make face out of setFaces
366  indirectPrimitivePatch bigFace
367  (
369  (
370  mesh_.faces(),
371  iter()
372  ),
373  mesh_.points()
374  );
375 
376  // Only store if -only one outside loop -which forms convex face
377  if (validFace(minConcaveCos, bigFace))
378  {
379  allFaceSets.append(iter());
380  }
381  }
382  }
383  }
384 
385  return allFaceSets.shrink();
386 }
387 
388 
390 (
391  const scalar featureCos,
392  const scalar minConcaveCos,
393  const labelHashSet& patchIDs
394 ) const
395 {
396  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
397 
398  // Pick up all cells on boundary
399  labelHashSet boundaryCells(mesh_.nFaces()-mesh_.nInternalFaces());
400 
401  forAllConstIter(labelHashSet, patchIDs, iter)
402  {
403  const polyPatch& patch = patches[iter.key()];
404 
405  if (!patch.coupled())
406  {
407  forAll(patch, i)
408  {
409  boundaryCells.insert(mesh_.faceOwner()[patch.start()+i]);
410  }
411  }
412  }
413 
414  return getMergeSets(featureCos, minConcaveCos, patchIDs, boundaryCells);
415 }
416 
417 
418 
420 (
421  const scalar featureCos,
422  const scalar minConcaveCos
423 ) const
424 {
425  const labelHashSet patchIDs(identity(mesh_.boundaryMesh().size()));
426 
427  return getMergeSets(featureCos, minConcaveCos, patchIDs);
428 }
429 
430 
431 // Gets outside edgeloop as a face
432 // - in same order as faces
433 // - in mesh vertex labels
435 (
436  const indirectPrimitivePatch& fp
437 )
438 {
439  if (fp.edgeLoops().size() != 1)
440  {
442  << "Multiple outside loops:" << fp.edgeLoops()
443  << abort(FatalError);
444  }
445 
446  // Get first boundary edge. Since guaranteed one edgeLoop when in here this
447  // edge must be on it.
448  label bEdgeI = fp.nInternalEdges();
449 
450  const edge& e = fp.edges()[bEdgeI];
451 
452  const labelList& eFaces = fp.edgeFaces()[bEdgeI];
453 
454  if (eFaces.size() != 1)
455  {
457  << "boundary edge:" << bEdgeI
458  << " points:" << fp.meshPoints()[e[0]]
459  << ' ' << fp.meshPoints()[e[1]]
460  << " on indirectPrimitivePatch has " << eFaces.size()
461  << " faces using it" << abort(FatalError);
462  }
463 
464 
465  // Outside loop
466  const labelList& outsideLoop = fp.edgeLoops()[0];
467 
468 
469  // Get order of edge e in outside loop
470  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
471 
472  // edgeLoopConsistent : true = edge in same order as outsideloop
473  // false = opposite order
474  bool edgeLoopConsistent = false;
475 
476  {
477  label index0 = findIndex(outsideLoop, e[0]);
478  label index1 = findIndex(outsideLoop, e[1]);
479 
480  if (index0 == -1 || index1 == -1)
481  {
483  << "Cannot find boundary edge:" << e
484  << " points:" << fp.meshPoints()[e[0]]
485  << ' ' << fp.meshPoints()[e[1]]
486  << " in edgeLoop:" << outsideLoop << abort(FatalError);
487  }
488  else if (index1 == outsideLoop.fcIndex(index0))
489  {
490  edgeLoopConsistent = true;
491  }
492  else if (index0 == outsideLoop.fcIndex(index1))
493  {
494  edgeLoopConsistent = false;
495  }
496  else
497  {
499  << "Cannot find boundary edge:" << e
500  << " points:" << fp.meshPoints()[e[0]]
501  << ' ' << fp.meshPoints()[e[1]]
502  << " on consecutive points in edgeLoop:"
503  << outsideLoop << abort(FatalError);
504  }
505  }
506 
507 
508  // Get face in local vertices
509  const face& localF = fp.localFaces()[eFaces[0]];
510 
511  // Get order of edge in localF
512  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
513 
514  // faceEdgeConsistent : true = edge in same order as localF
515  // false = opposite order
516  bool faceEdgeConsistent = false;
517 
518  {
519  // Find edge in face.
520  label index = findIndex(fp.faceEdges()[eFaces[0]], bEdgeI);
521 
522  if (index == -1)
523  {
525  << "Cannot find boundary edge:" << e
526  << " points:" << fp.meshPoints()[e[0]]
527  << ' ' << fp.meshPoints()[e[1]]
528  << " in face:" << eFaces[0]
529  << " edges:" << fp.faceEdges()[eFaces[0]]
530  << abort(FatalError);
531  }
532  else if (localF[index] == e[0] && localF.nextLabel(index) == e[1])
533  {
534  faceEdgeConsistent = true;
535  }
536  else if (localF[index] == e[1] && localF.nextLabel(index) == e[0])
537  {
538  faceEdgeConsistent = false;
539  }
540  else
541  {
543  << "Cannot find boundary edge:" << e
544  << " points:" << fp.meshPoints()[e[0]]
545  << ' ' << fp.meshPoints()[e[1]]
546  << " in face:" << eFaces[0] << " verts:" << localF
547  << abort(FatalError);
548  }
549  }
550 
551  // Get face in mesh points.
552  face meshFace(renumber(fp.meshPoints(), outsideLoop));
553 
554  if (faceEdgeConsistent != edgeLoopConsistent)
555  {
556  reverse(meshFace);
557  }
558  return meshFace;
559 }
560 
561 
563 (
564  const labelListList& faceSets,
565  polyTopoChange& meshMod
566 )
567 {
568  if (undoable_)
569  {
570  masterFace_.setSize(faceSets.size());
571  faceSetsVertices_.setSize(faceSets.size());
572  forAll(faceSets, setI)
573  {
574  const labelList& setFaces = faceSets[setI];
575 
576  masterFace_[setI] = setFaces[0];
577  faceSetsVertices_[setI] = UIndirectList<face>
578  (
579  mesh_.faces(),
580  setFaces
581  );
582  }
583  }
584 
585  // Running count of number of faces using a point
586  labelList nPointFaces(mesh_.nPoints(), 0);
587 
588  const labelListList& pointFaces = mesh_.pointFaces();
589 
590  forAll(pointFaces, pointi)
591  {
592  nPointFaces[pointi] = pointFaces[pointi].size();
593  }
594 
595  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
596 
597  forAll(faceSets, setI)
598  {
599  const labelList& setFaces = faceSets[setI];
600 
601  // Check
602  if (debug)
603  {
604  forAll(setFaces, i)
605  {
606  label patchi = patches.whichPatch(setFaces[i]);
607 
608  if (patchi == -1 || patches[patchi].coupled())
609  {
611  << "Can only merge non-coupled boundary faces"
612  << " but found internal or coupled face:"
613  << setFaces[i] << " in set " << setI
614  << abort(FatalError);
615  }
616  }
617  }
618 
619  // Make face out of setFaces
620  indirectPrimitivePatch bigFace
621  (
623  (
624  mesh_.faces(),
625  setFaces
626  ),
627  mesh_.points()
628  );
629 
630  // Get outside vertices (in local vertex numbering)
631  const labelListList& edgeLoops = bigFace.edgeLoops();
632 
633  if (edgeLoops.size() != 1)
634  {
636  << "Faces to-be-merged " << setFaces
637  << " do not form a single big face." << nl
638  << abort(FatalError);
639  }
640 
641 
642  // Delete all faces except master, whose face gets modified.
643 
644  // Modify master face
645  // ~~~~~~~~~~~~~~~~~~
646 
647  label masterFacei = setFaces[0];
648 
649  // Get outside face in mesh vertex labels
650  face outsideFace(getOutsideFace(bigFace));
651 
652  label zoneID = mesh_.faceZones().whichZone(masterFacei);
653 
654  bool zoneFlip = false;
655 
656  if (zoneID >= 0)
657  {
658  const faceZone& fZone = mesh_.faceZones()[zoneID];
659 
660  zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
661  }
662 
663  label patchi = mesh_.boundaryMesh().whichPatch(masterFacei);
664 
665  meshMod.setAction
666  (
668  (
669  outsideFace, // modified face
670  masterFacei, // label of face being modified
671  mesh_.faceOwner()[masterFacei], // owner
672  -1, // neighbour
673  false, // face flip
674  patchi, // patch for face
675  false, // remove from zone
676  zoneID, // zone for face
677  zoneFlip // face flip in zone
678  )
679  );
680 
681 
682  // Delete all non-master faces
683  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
684 
685  for (label i = 1; i < setFaces.size(); i++)
686  {
687  meshMod.setAction(polyRemoveFace(setFaces[i]));
688  }
689 
690 
691  // Mark unused points for deletion
692  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
693 
694  // Remove count of all faces combined
695  forAll(setFaces, i)
696  {
697  const face& f = mesh_.faces()[setFaces[i]];
698 
699  forAll(f, fp)
700  {
701  nPointFaces[f[fp]]--;
702  }
703  }
704  // But keep count on new outside face
705  forAll(outsideFace, fp)
706  {
707  nPointFaces[outsideFace[fp]]++;
708  }
709  }
710 
711 
712  // If one or more processors use the point it needs to be kept.
714  (
715  mesh_,
716  nPointFaces,
717  plusEqOp<label>(),
718  label(0) // null value
719  );
720 
721  // Remove all unused points. Store position if undoable.
722  if (!undoable_)
723  {
724  forAll(nPointFaces, pointi)
725  {
726  if (nPointFaces[pointi] == 0)
727  {
728  meshMod.setAction(polyRemovePoint(pointi));
729  }
730  }
731  }
732  else
733  {
734  // Count removed points
735  label n = 0;
736  forAll(nPointFaces, pointi)
737  {
738  if (nPointFaces[pointi] == 0)
739  {
740  n++;
741  }
742  }
743 
744  savedPointLabels_.setSize(n);
745  savedPoints_.setSize(n);
746  Map<label> meshToSaved(2*n);
747 
748  // Remove points and store position
749  n = 0;
750  forAll(nPointFaces, pointi)
751  {
752  if (nPointFaces[pointi] == 0)
753  {
754  meshMod.setAction(polyRemovePoint(pointi));
755 
756  savedPointLabels_[n] = pointi;
757  savedPoints_[n] = mesh_.points()[pointi];
758 
759  meshToSaved.insert(pointi, n);
760  n++;
761  }
762  }
763 
764  // Update stored vertex labels. Negative indices index into local points
765  forAll(faceSetsVertices_, setI)
766  {
767  faceList& setFaces = faceSetsVertices_[setI];
768 
769  forAll(setFaces, i)
770  {
771  face& f = setFaces[i];
772 
773  forAll(f, fp)
774  {
775  label pointi = f[fp];
776 
777  if (nPointFaces[pointi] == 0)
778  {
779  f[fp] = -meshToSaved[pointi]-1;
780  }
781  }
782  }
783  }
784  }
785 }
786 
787 
789 {
790  if (undoable_)
791  {
792  // Master face just renumbering of point labels
793  inplaceRenumber(map.reverseFaceMap(), masterFace_);
794 
795  // Stored faces refer to backed-up vertices (not changed)
796  // and normal vertices (need to be renumbered)
797  forAll(faceSetsVertices_, setI)
798  {
799  faceList& faces = faceSetsVertices_[setI];
800 
801  forAll(faces, i)
802  {
803  // Note: faces[i] can have negative elements.
804  face& f = faces[i];
805 
806  forAll(f, fp)
807  {
808  label pointi = f[fp];
809 
810  if (pointi >= 0)
811  {
812  f[fp] = map.reversePointMap()[pointi];
813 
814  if (f[fp] < 0)
815  {
817  << "In set " << setI << " at position " << i
818  << " with master face "
819  << masterFace_[setI] << nl
820  << "the points of the slave face " << faces[i]
821  << " don't exist anymore."
822  << abort(FatalError);
823  }
824  }
825  }
826  }
827  }
828  }
829 }
830 
831 
832 // Note that faces on coupled patches are never combined (or at least
833 // getMergeSets never picks them up. Thus the points on them will never get
834 // deleted since still used by the coupled faces. So never a risk of undoing
835 // things (faces/points) on coupled patches.
837 (
838  const labelList& masterFaces,
839  polyTopoChange& meshMod,
840  Map<label>& restoredPoints,
841  Map<label>& restoredFaces,
842  Map<label>& restoredCells
843 )
844 {
845  if (!undoable_)
846  {
848  << "Can only call setUnrefinement if constructed with"
849  << " unrefinement capability." << exit(FatalError);
850  }
851 
852 
853  // Restored points
854  labelList addedPoints(savedPoints_.size(), -1);
855 
856  // Invert set-to-master-face
857  Map<label> masterToSet(masterFace_.size());
858 
859  forAll(masterFace_, setI)
860  {
861  if (masterFace_[setI] >= 0)
862  {
863  masterToSet.insert(masterFace_[setI], setI);
864  }
865  }
866 
867  forAll(masterFaces, i)
868  {
869  label masterFacei = masterFaces[i];
870 
871  Map<label>::const_iterator iter = masterToSet.find(masterFacei);
872 
873  if (iter == masterToSet.end())
874  {
876  << "Master face " << masterFacei
877  << " is not the master of one of the merge sets"
878  << " or has already been merged"
879  << abort(FatalError);
880  }
881 
882  label setI = iter();
883 
884 
885  // Update faces of the merge setI for reintroduced vertices
886  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
887 
888  faceList& faces = faceSetsVertices_[setI];
889 
890  if (faces.empty())
891  {
893  << "Set " << setI << " with master face " << masterFacei
894  << " has already been merged." << abort(FatalError);
895  }
896 
897  forAll(faces, i)
898  {
899  face& f = faces[i];
900 
901  forAll(f, fp)
902  {
903  label pointi = f[fp];
904 
905  if (pointi < 0)
906  {
907  label localI = -pointi-1;
908 
909  if (addedPoints[localI] == -1)
910  {
911  // First occurrence of saved point. Reintroduce point
912  addedPoints[localI] = meshMod.setAction
913  (
915  (
916  savedPoints_[localI], // point
917  -1, // master point
918  -1, // zone for point
919  true // supports a cell
920  )
921  );
922  restoredPoints.insert
923  (
924  addedPoints[localI], // current point label
925  savedPointLabels_[localI] // point label when it
926  // was stored
927  );
928  }
929  f[fp] = addedPoints[localI];
930  }
931  }
932  }
933 
934 
935  // Restore
936  // ~~~~~~~
937 
938  label own = mesh_.faceOwner()[masterFacei];
939  label zoneID = mesh_.faceZones().whichZone(masterFacei);
940  bool zoneFlip = false;
941  if (zoneID >= 0)
942  {
943  const faceZone& fZone = mesh_.faceZones()[zoneID];
944  zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
945  }
946  label patchi = mesh_.boundaryMesh().whichPatch(masterFacei);
947 
948  if (mesh_.boundaryMesh()[patchi].coupled())
949  {
951  << "Master face " << masterFacei << " is on coupled patch "
952  << mesh_.boundaryMesh()[patchi].name()
953  << abort(FatalError);
954  }
955 
956  // Pout<< "Restoring new master face " << masterFacei
957  // << " to vertices " << faces[0] << endl;
958 
959  // Modify the master face.
960  meshMod.setAction
961  (
963  (
964  faces[0], // original face
965  masterFacei, // label of face
966  own, // owner
967  -1, // neighbour
968  false, // face flip
969  patchi, // patch for face
970  false, // remove from zone
971  zoneID, // zone for face
972  zoneFlip // face flip in zone
973  )
974  );
975  restoredFaces.insert(masterFacei, masterFacei);
976 
977  // Add the previously removed faces
978  for (label i = 1; i < faces.size(); i++)
979  {
980  // Pout<< "Restoring removed face with vertices " << faces[i]
981  // << endl;
982 
983  label facei = meshMod.setAction
984  (
986  (
987  faces[i], // vertices
988  own, // owner,
989  -1, // neighbour,
990  -1, // masterPointID,
991  -1, // masterEdgeID,
992  masterFacei, // masterFaceID,
993  false, // flipFaceFlux,
994  patchi, // patchID,
995  zoneID, // zoneID,
996  zoneFlip // zoneFlip
997  )
998  );
999  restoredFaces.insert(facei, masterFacei);
1000  }
1001 
1002  // Clear out restored set
1003  faceSetsVertices_[setI].clear();
1004  masterFace_[setI] = -1;
1005  }
1006 }
1007 
1008 
1009 // ************************************************************************* //
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:434
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:458
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:837
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:112
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:459
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:248
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.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
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
combineFaces(const polyMesh &mesh, const bool undoable=false)
Construct from mesh.
Definition: combineFaces.C:288
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:108
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.
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:305
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:208
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:435
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:313
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.
void setRefinement(const labelListList &, polyTopoChange &)
Play commands into polyTopoChange to combine faces. Gets.
Definition: combineFaces.C:563
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
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
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:303
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:788
A List with indirect addressing.
Definition: IndirectList.H:101
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:490
Namespace for OpenFOAM.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.