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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "combineFaces.H"
30 #include "polyMesh.H"
31 #include "polyTopoChange.H"
32 #include "polyRemoveFace.H"
33 #include "polyAddFace.H"
34 #include "polyModifyFace.H"
35 #include "polyRemovePoint.H"
36 #include "polyAddPoint.H"
37 #include "syncTools.H"
38 #include "meshTools.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(combineFaces, 0);
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 // Test face for (almost) convexeness. Allows certain convexity before
51 // complaining.
52 // For every two consecutive edges calculate the normal. If it differs too
53 // much from the average face normal complain.
54 bool Foam::combineFaces::convexFace
55 (
56  const scalar minConcaveCos,
57  const pointField& points,
58  const face& f
59 )
60 {
61  // Get outwards pointing normal of f, only the sign matters.
62  const vector areaNorm = f.areaNormal(points);
63 
64  // Get edge from f[0] to f[size-1];
65  vector ePrev(points[f.first()] - points[f.last()]);
66  scalar magEPrev = mag(ePrev);
67  ePrev /= magEPrev + VSMALL;
68 
69  forAll(f, fp0)
70  {
71  // Get vertex after fp
72  label fp1 = f.fcIndex(fp0);
73 
74  // Normalized vector between two consecutive points
75  vector e10(points[f[fp1]] - points[f[fp0]]);
76  scalar magE10 = mag(e10);
77  e10 /= magE10 + VSMALL;
78 
79  if (magEPrev > SMALL && magE10 > SMALL)
80  {
81  vector edgeNormal = ePrev ^ e10;
82 
83  if ((edgeNormal & areaNorm) < 0)
84  {
85  // Concave. Check angle.
86  if ((ePrev & e10) < minConcaveCos)
87  {
88  return false;
89  }
90  }
91  }
92 
93  ePrev = e10;
94  magEPrev = magE10;
95  }
96 
97  // Not a single internal angle is concave so face is convex.
98  return true;
99 }
100 
101 
102 // Determines if set of faces is valid to collapse into single face.
103 bool Foam::combineFaces::validFace
104 (
105  const scalar minConcaveCos,
106  const indirectPrimitivePatch& bigFace
107 )
108 {
109  // Get outside vertices (in local vertex numbering)
110  const labelListList& edgeLoops = bigFace.edgeLoops();
111 
112  if (edgeLoops.size() > 1)
113  {
114  return false;
115  }
116 
117  bool isNonManifold = bigFace.checkPointManifold(false, nullptr);
118  if (isNonManifold)
119  {
120  return false;
121  }
122 
123  // Check for convexness
124  face f(getOutsideFace(bigFace));
125 
126  return convexFace(minConcaveCos, bigFace.points(), f);
127 }
128 
129 
130 void Foam::combineFaces::regioniseFaces
131 (
132  const scalar minCos,
133  const bool mergeAcrossPatches,
134  const label celli,
135  const labelList& cEdges,
136  Map<label>& faceRegion
137 ) const
138 {
139  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
140 
141  forAll(cEdges, i)
142  {
143  const label edgeI = cEdges[i];
144 
145  label f0, f1;
146  meshTools::getEdgeFaces(mesh_, celli, edgeI, f0, f1);
147 
148  const vector& a0 = mesh_.faceAreas()[f0];
149  const vector& a1 = mesh_.faceAreas()[f1];
150 
151  const label p0 = patches.whichPatch(f0);
152  const label p1 = patches.whichPatch(f1);
153 
154  // Face can be merged if
155  // - small angle
156  // - mergeAcrossPatches=false : same non-constraint patch
157  // - mergeAcrossPatches=true : always (if non-constraint patch)
158  // (this logic could be extended to e.g. merge faces on symm plane
159  // if they have similar normals. But there might be lots of other
160  // constraints which disallow merging so this decision ideally should
161  // be up to patch type)
162  if
163  (
164  p0 != -1
165  && p1 != -1
166  && !(
169  )
170  )
171  {
172  if (!mergeAcrossPatches && (p0 != p1))
173  {
174  continue;
175  }
176 
177  const vector f0Normal = normalised(a0);
178  const vector f1Normal = normalised(a1);
179 
180  if ((f0Normal & f1Normal) > minCos)
181  {
182  const label region0 = faceRegion.lookup(f0, -1);
183  const label region1 = faceRegion.lookup(f1, -1);
184 
185  if (region0 == -1)
186  {
187  if (region1 == -1)
188  {
189  const label useRegion = faceRegion.size();
190  faceRegion.insert(f0, useRegion);
191  faceRegion.insert(f1, useRegion);
192  }
193  else
194  {
195  faceRegion.insert(f0, region1);
196  }
197  }
198  else
199  {
200  if (region1 == -1)
201  {
202  faceRegion.insert(f1, region0);
203  }
204  else if (region0 != region1)
205  {
206  // Merge the two regions
207  const label useRegion = min(region0, region1);
208  const label freeRegion = max(region0, region1);
209 
210  forAllIters(faceRegion, iter)
211  {
212  if (iter.val() == freeRegion)
213  {
214  iter.val() = useRegion;
215  }
216  }
217  }
218  }
219  }
220  }
221  }
222 }
223 
224 
225 bool Foam::combineFaces::faceNeighboursValid
226 (
227  const label celli,
228  const Map<label>& faceRegion
229 ) const
230 {
231  if (faceRegion.size() <= 1)
232  {
233  return true;
234  }
235 
236  const cell& cFaces = mesh_.cells()[celli];
237 
238  DynamicList<label> storage;
239 
240  // Test for face collapsing to edge since too many neighbours merged.
241  forAll(cFaces, cFacei)
242  {
243  label facei = cFaces[cFacei];
244 
245  if (!faceRegion.found(facei))
246  {
247  const labelList& fEdges = mesh_.faceEdges(facei, storage);
248 
249  // Count number of remaining faces neighbouring facei. This has
250  // to be 3 or more.
251 
252  // Unregioned neighbouring faces
253  DynamicList<label> neighbourFaces(cFaces.size());
254  // Regioned neighbouring faces
255  labelHashSet neighbourRegions;
256 
257  forAll(fEdges, i)
258  {
259  const label edgeI = fEdges[i];
260  label nbrI = meshTools::otherFace(mesh_, celli, facei, edgeI);
261 
262  const auto iter = faceRegion.cfind(nbrI);
263 
264  if (iter.found())
265  {
266  neighbourRegions.insert(iter.val());
267  }
268  else
269  {
270  neighbourFaces.appendUniq(nbrI);
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
287 Foam::combineFaces::combineFaces
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& boundaryCells,
309  const bool mergeAcrossPatches
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.
316  DynamicList<label> storage;
317 
318  // On all cells regionise the faces
319  for (const label celli : boundaryCells)
320  {
321  const cell& cFaces = mesh_.cells()[celli];
322 
323  const labelList& cEdges = mesh_.cellEdges(celli, set, storage);
324 
325  // Region per face
326  Map<label> faceRegion(cFaces.size());
327  regioniseFaces
328  (
329  featureCos,
330  mergeAcrossPatches,
331  celli,
332  cEdges,
333  faceRegion
334  );
335 
336  // Now we have in faceRegion for every face the region with planar
337  // face sharing the same region. We now check whether the resulting
338  // sets cause a face
339  // - to become a set of edges since too many faces are merged.
340  // - to become convex
341 
342  if (faceNeighboursValid(celli, faceRegion))
343  {
344  // Create region-to-faces addressing
345  Map<labelList> regionToFaces(faceRegion.size());
346 
347  forAllConstIters(faceRegion, iter)
348  {
349  const label facei = iter.key();
350  const label region = iter.val();
351 
352  auto regionFnd = regionToFaces.find(region);
353 
354  if (regionFnd.found())
355  {
356  labelList& setFaces = regionFnd();
357  label sz = setFaces.size();
358  setFaces.setSize(sz+1);
359  setFaces[sz] = facei;
360  }
361  else
362  {
363  regionToFaces.insert(region, labelList(1, facei));
364  }
365  }
366 
367  // For every set check if it forms a valid convex face
368 
369  forAllIters(regionToFaces, iter)
370  {
371  // Make face out of setFaces
372  indirectPrimitivePatch bigFace
373  (
375  (
376  mesh_.faces(),
377  iter.val()
378  ),
379  mesh_.points()
380  );
381 
382  // Only store if -only one outside loop -which forms convex face
383  if (validFace(minConcaveCos, bigFace))
384  {
385  labelList& faceIDs = iter.val();
386 
387  // For cross-patch merging we want to make the
388  // largest face the one to decide the final patch
389  // (i.e. master face)
390  if (mergeAcrossPatches)
391  {
392  const vectorField& areas = mesh_.faceAreas();
393 
394  label maxIndex = 0;
395  scalar maxMagSqr = magSqr(areas[faceIDs[0]]);
396  for (label i = 1; i < faceIDs.size(); ++i)
397  {
398  const scalar a2 = magSqr(areas[faceIDs[i]]);
399  if (a2 > maxMagSqr)
400  {
401  maxMagSqr = a2;
402  maxIndex = i;
403  }
404  }
405  if (maxIndex != 0)
406  {
407  std::swap(faceIDs[0], faceIDs[maxIndex]);
408  }
409  }
410 
411  allFaceSets.append(faceIDs);
412  }
413  }
414  }
415  }
416 
417  return allFaceSets.shrink();
418 }
419 
420 
422 (
423  const scalar featureCos,
424  const scalar minConcaveCos,
425  const bool mergeAcrossPatches
426 ) const
427 {
428  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
429 
430  // Pick up all cells on boundary
431  labelHashSet boundaryCells(mesh_.nBoundaryFaces());
432 
433  forAll(patches, patchi)
434  {
435  const polyPatch& patch = patches[patchi];
436 
437  if (!patch.coupled())
438  {
439  forAll(patch, i)
440  {
441  boundaryCells.insert(mesh_.faceOwner()[patch.start()+i]);
442  }
443  }
444  }
445 
446  return getMergeSets
447  (
448  featureCos,
449  minConcaveCos,
450  boundaryCells,
451  mergeAcrossPatches
452  );
453 }
454 
455 
456 // Gets outside edgeloop as a face
457 // - in same order as faces
458 // - in mesh vertex labels
460 (
461  const indirectPrimitivePatch& fp
462 )
463 {
464  if (fp.edgeLoops().size() != 1)
465  {
467  << "Multiple outside loops:" << fp.edgeLoops()
468  << abort(FatalError);
469  }
470 
471  // Get first boundary edge. Since guaranteed one edgeLoop when in here this
472  // edge must be on it.
473  label bEdgeI = fp.nInternalEdges();
474 
475  const edge& e = fp.edges()[bEdgeI];
476 
477  const labelList& eFaces = fp.edgeFaces()[bEdgeI];
478 
479  if (eFaces.size() != 1)
480  {
482  << "boundary edge:" << bEdgeI
483  << " points:" << fp.meshPoints()[e[0]]
484  << ' ' << fp.meshPoints()[e[1]]
485  << " on indirectPrimitivePatch has " << eFaces.size()
486  << " faces using it" << abort(FatalError);
487  }
488 
489 
490  // Outside loop
491  const labelList& outsideLoop = fp.edgeLoops()[0];
492 
493 
494  // Get order of edge e in outside loop
495  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
496 
497  // edgeLoopConsistent : true = edge in same order as outsideloop
498  // false = opposite order
499  bool edgeLoopConsistent = false;
500 
501  {
502  label index0 = outsideLoop.find(e[0]);
503  label index1 = outsideLoop.find(e[1]);
504 
505  if (index0 == -1 || index1 == -1)
506  {
508  << "Cannot find boundary edge:" << e
509  << " points:" << fp.meshPoints()[e[0]]
510  << ' ' << fp.meshPoints()[e[1]]
511  << " in edgeLoop:" << outsideLoop << abort(FatalError);
512  }
513  else if (index1 == outsideLoop.fcIndex(index0))
514  {
515  edgeLoopConsistent = true;
516  }
517  else if (index0 == outsideLoop.fcIndex(index1))
518  {
519  edgeLoopConsistent = false;
520  }
521  else
522  {
524  << "Cannot find boundary edge:" << e
525  << " points:" << fp.meshPoints()[e[0]]
526  << ' ' << fp.meshPoints()[e[1]]
527  << " on consecutive points in edgeLoop:"
528  << outsideLoop << abort(FatalError);
529  }
530  }
531 
532 
533  // Get face in local vertices
534  const face& localF = fp.localFaces()[eFaces[0]];
535 
536  // Get order of edge in localF
537  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
538 
539  // faceEdgeConsistent : true = edge in same order as localF
540  // false = opposite order
541  bool faceEdgeConsistent = false;
542 
543  {
544  // Find edge in face.
545  label index = fp.faceEdges()[eFaces[0]].find(bEdgeI);
546 
547  if (index == -1)
548  {
550  << "Cannot find boundary edge:" << e
551  << " points:" << fp.meshPoints()[e[0]]
552  << ' ' << fp.meshPoints()[e[1]]
553  << " in face:" << eFaces[0]
554  << " edges:" << fp.faceEdges()[eFaces[0]]
555  << abort(FatalError);
556  }
557  else if (localF[index] == e[0] && localF.nextLabel(index) == e[1])
558  {
559  faceEdgeConsistent = true;
560  }
561  else if (localF[index] == e[1] && localF.nextLabel(index) == e[0])
562  {
563  faceEdgeConsistent = false;
564  }
565  else
566  {
568  << "Cannot find boundary edge:" << e
569  << " points:" << fp.meshPoints()[e[0]]
570  << ' ' << fp.meshPoints()[e[1]]
571  << " in face:" << eFaces[0] << " verts:" << localF
572  << abort(FatalError);
573  }
574  }
575 
576  // Get face in mesh points.
577  face meshFace(renumber(fp.meshPoints(), outsideLoop));
578 
579  if (faceEdgeConsistent != edgeLoopConsistent)
580  {
581  reverse(meshFace);
582  }
583  return meshFace;
584 }
585 
586 
588 (
589  const labelListList& faceSets,
590  polyTopoChange& meshMod
591 )
592 {
593  if (undoable_)
594  {
595  masterFace_.setSize(faceSets.size());
596  faceSetsVertices_.setSize(faceSets.size());
597  forAll(faceSets, setI)
598  {
599  const labelList& setFaces = faceSets[setI];
600 
601  masterFace_[setI] = setFaces[0];
602  faceSetsVertices_[setI] = UIndirectList<face>
603  (
604  mesh_.faces(),
605  setFaces
606  );
607  }
608  }
609 
610  // Running count of number of faces using a point
611  labelList nPointFaces(mesh_.nPoints(), Zero);
612 
613  const labelListList& pointFaces = mesh_.pointFaces();
614 
615  forAll(pointFaces, pointi)
616  {
617  nPointFaces[pointi] = pointFaces[pointi].size();
618  }
619 
620  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
621 
622  forAll(faceSets, setI)
623  {
624  const labelList& setFaces = faceSets[setI];
625 
626  // Check
627  if (debug)
628  {
629  forAll(setFaces, i)
630  {
631  label patchi = patches.whichPatch(setFaces[i]);
632 
633  if (patchi == -1 || patches[patchi].coupled())
634  {
636  << "Can only merge non-coupled boundary faces"
637  << " but found internal or coupled face:"
638  << setFaces[i] << " in set " << setI
639  << abort(FatalError);
640  }
641  }
642  }
643 
644  // Make face out of setFaces
645  indirectPrimitivePatch bigFace
646  (
648  (
649  mesh_.faces(),
650  setFaces
651  ),
652  mesh_.points()
653  );
654 
655  // Get outside vertices (in local vertex numbering)
656  const labelListList& edgeLoops = bigFace.edgeLoops();
657 
658  if (edgeLoops.size() != 1)
659  {
661  << "Faces to-be-merged " << setFaces
662  << " do not form a single big face." << nl
663  << abort(FatalError);
664  }
665 
666 
667  // Delete all faces except master, whose face gets modified.
668 
669  // Modify master face
670  // ~~~~~~~~~~~~~~~~~~
671 
672  label masterFacei = setFaces[0];
673 
674  // Get outside face in mesh vertex labels
675  face outsideFace(getOutsideFace(bigFace));
676 
677  label zoneID = mesh_.faceZones().whichZone(masterFacei);
678 
679  bool zoneFlip = false;
680 
681  if (zoneID >= 0)
682  {
683  const faceZone& fZone = mesh_.faceZones()[zoneID];
684 
685  zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
686  }
687 
688  label patchi = mesh_.boundaryMesh().whichPatch(masterFacei);
689 
690  meshMod.setAction
691  (
693  (
694  outsideFace, // modified face
695  masterFacei, // label of face being modified
696  mesh_.faceOwner()[masterFacei], // owner
697  -1, // neighbour
698  false, // face flip
699  patchi, // patch for face
700  false, // remove from zone
701  zoneID, // zone for face
702  zoneFlip // face flip in zone
703  )
704  );
705 
706 
707  // Delete all non-master faces
708  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
709 
710  for (label i = 1; i < setFaces.size(); i++)
711  {
712  meshMod.setAction(polyRemoveFace(setFaces[i]));
713  }
714 
715 
716  // Mark unused points for deletion
717  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
718 
719  // Remove count of all faces combined
720  forAll(setFaces, i)
721  {
722  const face& f = mesh_.faces()[setFaces[i]];
723 
724  forAll(f, fp)
725  {
726  nPointFaces[f[fp]]--;
727  }
728  }
729  // But keep count on new outside face
730  forAll(outsideFace, fp)
731  {
732  nPointFaces[outsideFace[fp]]++;
733  }
734  }
735 
736 
737  // If one or more processors use the point it needs to be kept.
739  (
740  mesh_,
741  nPointFaces,
742  plusEqOp<label>(),
743  label(0) // null value
744  );
745 
746  // Remove all unused points. Store position if undoable.
747  if (!undoable_)
748  {
749  forAll(nPointFaces, pointi)
750  {
751  if (nPointFaces[pointi] == 0)
752  {
753  meshMod.setAction(polyRemovePoint(pointi));
754  }
755  }
756  }
757  else
758  {
759  // Count removed points
760  label n = 0;
761  forAll(nPointFaces, pointi)
762  {
763  if (nPointFaces[pointi] == 0)
764  {
765  n++;
766  }
767  }
768 
769  savedPointLabels_.setSize(n);
770  savedPoints_.setSize(n);
771  Map<label> meshToSaved(2*n);
772 
773  // Remove points and store position
774  n = 0;
775  forAll(nPointFaces, pointi)
776  {
777  if (nPointFaces[pointi] == 0)
778  {
779  meshMod.setAction(polyRemovePoint(pointi));
780 
781  savedPointLabels_[n] = pointi;
782  savedPoints_[n] = mesh_.points()[pointi];
783 
784  meshToSaved.insert(pointi, n);
785  n++;
786  }
787  }
788 
789  // Update stored vertex labels. Negative indices index into local points
790  forAll(faceSetsVertices_, setI)
791  {
792  faceList& setFaces = faceSetsVertices_[setI];
793 
794  forAll(setFaces, i)
795  {
796  face& f = setFaces[i];
797 
798  forAll(f, fp)
799  {
800  label pointi = f[fp];
801 
802  if (nPointFaces[pointi] == 0)
803  {
804  f[fp] = -meshToSaved[pointi]-1;
805  }
806  }
807  }
808  }
809  }
810 }
811 
812 
814 {
815  if (undoable_)
816  {
817  // Master face just renumbering of point labels
818  inplaceRenumber(map.reverseFaceMap(), masterFace_);
819 
820  // Stored faces refer to backed-up vertices (not changed)
821  // and normal vertices (need to be renumbered)
822  forAll(faceSetsVertices_, setI)
823  {
824  faceList& faces = faceSetsVertices_[setI];
825 
826  forAll(faces, i)
827  {
828  // Note: faces[i] can have negative elements.
829  face& f = faces[i];
830 
831  forAll(f, fp)
832  {
833  label pointi = f[fp];
834 
835  if (pointi >= 0)
836  {
837  f[fp] = map.reversePointMap()[pointi];
838 
839  if (f[fp] < 0)
840  {
842  << "In set " << setI << " at position " << i
843  << " with master face "
844  << masterFace_[setI] << nl
845  << "the points of the slave face " << faces[i]
846  << " don't exist anymore."
847  << abort(FatalError);
848  }
849  }
850  }
851  }
852  }
853  }
854 }
855 
856 
857 // Note that faces on coupled patches are never combined (or at least
858 // getMergeSets never picks them up. Thus the points on them will never get
859 // deleted since still used by the coupled faces. So never a risk of undoing
860 // things (faces/points) on coupled patches.
862 (
863  const labelList& masterFaces,
864  polyTopoChange& meshMod,
865  Map<label>& restoredPoints,
866  Map<label>& restoredFaces,
867  Map<label>& restoredCells
868 )
869 {
870  if (!undoable_)
871  {
873  << "Can only call setUnrefinement if constructed with"
874  << " unrefinement capability." << exit(FatalError);
875  }
876 
877 
878  // Restored points
879  labelList addedPoints(savedPoints_.size(), -1);
880 
881  // Invert set-to-master-face
882  Map<label> masterToSet(masterFace_.size());
883 
884  forAll(masterFace_, setI)
885  {
886  if (masterFace_[setI] >= 0)
887  {
888  masterToSet.insert(masterFace_[setI], setI);
889  }
890  }
891 
892  forAll(masterFaces, i)
893  {
894  const label masterFacei = masterFaces[i];
895 
896  const auto iter = masterToSet.cfind(masterFacei);
897 
898  if (!iter.found())
899  {
901  << "Master face " << masterFacei
902  << " is not the master of one of the merge sets"
903  << " or has already been merged"
904  << abort(FatalError);
905  }
906 
907  const label setI = iter.val();
908 
909 
910  // Update faces of the merge setI for reintroduced vertices
911  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
912 
913  faceList& faces = faceSetsVertices_[setI];
914 
915  if (faces.empty())
916  {
918  << "Set " << setI << " with master face " << masterFacei
919  << " has already been merged." << abort(FatalError);
920  }
921 
922  forAll(faces, i)
923  {
924  face& f = faces[i];
925 
926  forAll(f, fp)
927  {
928  label pointi = f[fp];
929 
930  if (pointi < 0)
931  {
932  label localI = -pointi-1;
933 
934  if (addedPoints[localI] == -1)
935  {
936  // First occurrence of saved point. Reintroduce point
937  addedPoints[localI] = meshMod.setAction
938  (
940  (
941  savedPoints_[localI], // point
942  -1, // master point
943  -1, // zone for point
944  true // supports a cell
945  )
946  );
947  restoredPoints.insert
948  (
949  addedPoints[localI], // current point label
950  savedPointLabels_[localI] // point label when it
951  // was stored
952  );
953  }
954  f[fp] = addedPoints[localI];
955  }
956  }
957  }
958 
959 
960  // Restore
961  // ~~~~~~~
962 
963  label own = mesh_.faceOwner()[masterFacei];
964  label zoneID = mesh_.faceZones().whichZone(masterFacei);
965  bool zoneFlip = false;
966  if (zoneID >= 0)
967  {
968  const faceZone& fZone = mesh_.faceZones()[zoneID];
969  zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
970  }
971  label patchi = mesh_.boundaryMesh().whichPatch(masterFacei);
972 
973  if (mesh_.boundaryMesh()[patchi].coupled())
974  {
976  << "Master face " << masterFacei << " is on coupled patch "
977  << mesh_.boundaryMesh()[patchi].name()
978  << abort(FatalError);
979  }
980 
981  //Pout<< "Restoring new master face " << masterFacei
982  // << " to vertices " << faces[0] << endl;
983 
984  // Modify the master face.
985  meshMod.setAction
986  (
988  (
989  faces[0], // original face
990  masterFacei, // label of face
991  own, // owner
992  -1, // neighbour
993  false, // face flip
994  patchi, // patch for face
995  false, // remove from zone
996  zoneID, // zone for face
997  zoneFlip // face flip in zone
998  )
999  );
1000  restoredFaces.insert(masterFacei, masterFacei);
1001 
1002  // Add the previously removed faces
1003  for (label i = 1; i < faces.size(); i++)
1004  {
1005  //Pout<< "Restoring removed face with vertices " << faces[i]
1006  // << endl;
1007 
1008  label facei = meshMod.setAction
1009  (
1010  polyAddFace
1011  (
1012  faces[i], // vertices
1013  own, // owner,
1014  -1, // neighbour,
1015  -1, // masterPointID,
1016  -1, // masterEdgeID,
1017  masterFacei, // masterFaceID,
1018  false, // flipFaceFlux,
1019  patchi, // patchID,
1020  zoneID, // zoneID,
1021  zoneFlip // zoneFlip
1022  )
1023  );
1024  restoredFaces.insert(facei, masterFacei);
1025  }
1026 
1027  // Clear out restored set
1028  faceSetsVertices_[setI].clear();
1029  masterFace_[setI] = -1;
1030  }
1031 }
1032 
1033 
1034 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::reverse
void reverse(UList< T > &list, const label n)
Definition: UListI.H:449
meshTools.H
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::combineFaces::getOutsideFace
static face getOutsideFace(const indirectPrimitivePatch &)
Gets outside of patch as a face (in mesh point labels)
Definition: combineFaces.C:460
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatch.C:262
polyRemovePoint.H
Foam::faceZone::flipMap
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:272
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatch.C:183
Foam::polyRemovePoint
Class containing data for point removal.
Definition: polyRemovePoint.H:48
Foam::meshTools::otherFace
label otherFace(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:555
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::combineFaces::getMergeSets
labelListList getMergeSets(const scalar featureCos, const scalar minConcaveCos, const labelHashSet &boundaryCells, const bool mergeAcrossPatches=false) const
Extract lists of all (non-coupled) boundary faces on selected.
Definition: combineFaces.C:305
polyAddFace.H
polyRemoveFace.H
polyTopoChange.H
combineFaces.H
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:99
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::Map< label >
Foam::polyModifyFace
Class describing modification of a face.
Definition: polyModifyFace.H:50
Foam::PrimitivePatch::faceEdges
const labelListList & faceEdges() const
Return face-edge addressing.
Definition: PrimitivePatch.C:275
polyMesh.H
Foam::HashSet< label, Hash< label > >
syncTools.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:61
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:721
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::maxMagSqr
Type maxMagSqr(const UList< Type > &f)
Definition: FieldFunctions.C:417
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::Field< vector >
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::combineFaces::updateMesh
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: combineFaces.C:813
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::combineFaces::setUnrefinement
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:862
Foam::meshTools::getEdgeFaces
void getEdgeFaces(const primitiveMesh &mesh, 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:479
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:56
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:812
polyAddPoint.H
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::polyPatch::constraintType
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:277
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::PrimitivePatch::localFaces
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatch.C:317
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:372
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
polyModifyFace.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:501
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::indirectPrimitivePatch
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
Definition: indirectPrimitivePatch.H:49
Foam::PrimitivePatch::nInternalEdges
label nInternalEdges() const
Number of internal edges.
Definition: PrimitivePatch.C:214
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::renumber
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:37
Foam::constant::atomic::a0
const dimensionedScalar a0
Bohr radius: default SI units: [m].
Foam::face::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:175
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::polyTopoChange::setAction
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
Definition: polyTopoChange.C:2481
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mapPolyMesh::reversePointMap
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:469
Foam::polyRemoveFace
Class containing data for face removal.
Definition: polyRemoveFace.H:48
Foam::PrimitivePatch::edgeLoops
const labelListList & edgeLoops() const
Return list of closed loops of boundary vertices.
Definition: PrimitivePatchEdgeLoops.C:146
Foam::List< labelList >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::polyAddPoint
Class containing data for point addition.
Definition: polyAddPoint.H:49
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::polyAddFace
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:53
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::combineFaces::setRefinement
void setRefinement(const labelListList &, polyTopoChange &)
Play commands into polyTopoChange to combine faces. Gets.
Definition: combineFaces.C:588
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::plusEqOp
Definition: ops.H:72
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:330
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79