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