removeFaces.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 "removeFaces.H"
27 #include "polyMesh.H"
28 #include "polyTopoChange.H"
29 #include "meshTools.H"
30 #include "polyModifyFace.H"
31 #include "polyRemoveFace.H"
32 #include "polyRemoveCell.H"
33 #include "polyRemovePoint.H"
34 #include "syncTools.H"
35 #include "OFstream.H"
36 #include "indirectPrimitivePatch.H"
37 #include "Time.H"
38 #include "faceSet.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 defineTypeNameAndDebug(removeFaces, 0);
46 
47 }
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 // Changes region of connected set of cells. Can be recursive since hopefully
52 // only small area of faces removed in one go.
53 void Foam::removeFaces::changeCellRegion
54 (
55  const label celli,
56  const label oldRegion,
57  const label newRegion,
58  labelList& cellRegion
59 ) const
60 {
61  if (cellRegion[celli] == oldRegion)
62  {
63  cellRegion[celli] = newRegion;
64 
65  // Step to neighbouring cells
66 
67  const labelList& cCells = mesh_.cellCells()[celli];
68 
69  forAll(cCells, i)
70  {
71  changeCellRegion(cCells[i], oldRegion, newRegion, cellRegion);
72  }
73  }
74 }
75 
76 
77 // Changes region of connected set of faces. Returns number of changed faces.
78 Foam::label Foam::removeFaces::changeFaceRegion
79 (
80  const labelList& cellRegion,
81  const boolList& removedFace,
82  const labelList& nFacesPerEdge,
83  const label facei,
84  const label newRegion,
85  const labelList& fEdges,
86  labelList& faceRegion
87 ) const
88 {
89  label nChanged = 0;
90 
91  if (faceRegion[facei] == -1 && !removedFace[facei])
92  {
93  faceRegion[facei] = newRegion;
94 
95  nChanged = 1;
96 
97  // Storage for on-the-fly addressing
98  DynamicList<label> fe;
99  DynamicList<label> ef;
100 
101  // Step to neighbouring faces across edges that will get removed
102  forAll(fEdges, i)
103  {
104  label edgeI = fEdges[i];
105 
106  if (nFacesPerEdge[edgeI] >= 0 && nFacesPerEdge[edgeI] <= 2)
107  {
108  const labelList& eFaces = mesh_.edgeFaces(edgeI, ef);
109 
110  forAll(eFaces, j)
111  {
112  label nbrFacei = eFaces[j];
113 
114  const labelList& fEdges1 = mesh_.faceEdges(nbrFacei, fe);
115 
116  nChanged += changeFaceRegion
117  (
118  cellRegion,
119  removedFace,
120  nFacesPerEdge,
121  nbrFacei,
122  newRegion,
123  fEdges1,
124  faceRegion
125  );
126  }
127  }
128  }
129  }
130  return nChanged;
131 }
132 
133 
134 // Mark all faces affected in any way by
135 // - removal of cells
136 // - removal of faces
137 // - removal of edges
138 // - removal of points
139 Foam::boolList Foam::removeFaces::getFacesAffected
140 (
141  const labelList& cellRegion,
142  const labelList& cellRegionMaster,
143  const labelList& facesToRemove,
144  const labelHashSet& edgesToRemove,
145  const labelHashSet& pointsToRemove
146 ) const
147 {
148  boolList affectedFace(mesh_.nFaces(), false);
149 
150  // Mark faces affected by removal of cells
151  forAll(cellRegion, celli)
152  {
153  label region = cellRegion[celli];
154 
155  if (region != -1 && (celli != cellRegionMaster[region]))
156  {
157  const labelList& cFaces = mesh_.cells()[celli];
158 
159  forAll(cFaces, cFacei)
160  {
161  affectedFace[cFaces[cFacei]] = true;
162  }
163  }
164  }
165 
166  // Mark faces affected by removal of face.
167  forAll(facesToRemove, i)
168  {
169  affectedFace[facesToRemove[i]] = true;
170  }
171 
172  // Mark faces affected by removal of edges
173  forAllConstIter(labelHashSet, edgesToRemove, iter)
174  {
175  const labelList& eFaces = mesh_.edgeFaces(iter.key());
176 
177  forAll(eFaces, eFacei)
178  {
179  affectedFace[eFaces[eFacei]] = true;
180  }
181  }
182 
183  // Mark faces affected by removal of points
184  forAllConstIter(labelHashSet, pointsToRemove, iter)
185  {
186  label pointi = iter.key();
187 
188  const labelList& pFaces = mesh_.pointFaces()[pointi];
189 
190  forAll(pFaces, pFacei)
191  {
192  affectedFace[pFaces[pFacei]] = true;
193  }
194  }
195  return affectedFace;
196 }
197 
198 
199 void Foam::removeFaces::writeOBJ
200 (
201  const indirectPrimitivePatch& fp,
202  const fileName& fName
203 )
204 {
205  OFstream str(fName);
206  Pout<< "removeFaces::writeOBJ : Writing faces to file "
207  << str.name() << endl;
208 
209  const pointField& localPoints = fp.localPoints();
210 
211  forAll(localPoints, i)
212  {
213  meshTools::writeOBJ(str, localPoints[i]);
214  }
215 
216  const faceList& localFaces = fp.localFaces();
217 
218  forAll(localFaces, i)
219  {
220  const face& f = localFaces[i];
221 
222  str<< 'f';
223 
224  forAll(f, fp)
225  {
226  str<< ' ' << f[fp]+1;
227  }
228  str<< nl;
229  }
230 }
231 
232 
233 // Inserts commands to merge faceLabels into one face.
234 void Foam::removeFaces::mergeFaces
235 (
236  const labelList& cellRegion,
237  const labelList& cellRegionMaster,
238  const labelHashSet& pointsToRemove,
239  const labelList& faceLabels,
240  polyTopoChange& meshMod
241 ) const
242 {
243  // Construct addressing engine from faceLabels (in order of faceLabels as
244  // well)
246  (
247  IndirectList<face>
248  (
249  mesh_.faces(),
250  faceLabels
251  ),
252  mesh_.points()
253  );
254 
255  // Get outside vertices (in local vertex numbering)
256 
257  if (fp.edgeLoops().size() != 1)
258  {
259  writeOBJ(fp, mesh_.time().path()/"facesToBeMerged.obj");
261  << "Cannot merge faces " << faceLabels
262  << " into single face since outside vertices " << fp.edgeLoops()
263  << " do not form single loop but form " << fp.edgeLoops().size()
264  << " loops instead." << abort(FatalError);
265  }
266 
267  const labelList& edgeLoop = fp.edgeLoops()[0];
268 
269  // Get outside vertices in order of one of the faces in faceLabels.
270  // (this becomes the master face)
271  // Find the first face that uses edgeLoop[0] and edgeLoop[1] as consecutive
272  // vertices.
273 
274  label masterIndex = -1;
275  bool reverseLoop = false;
276 
277  const labelList& pFaces = fp.pointFaces()[edgeLoop[0]];
278 
279  // Find face among pFaces which uses edgeLoop[1]
280  forAll(pFaces, i)
281  {
282  label facei = pFaces[i];
283 
284  const face& f = fp.localFaces()[facei];
285 
286  label index1 = findIndex(f, edgeLoop[1]);
287 
288  if (index1 != -1)
289  {
290  // Check whether consecutive to edgeLoop[0]
291  label index0 = findIndex(f, edgeLoop[0]);
292 
293  if (index0 != -1)
294  {
295  if (index1 == f.fcIndex(index0))
296  {
297  masterIndex = facei;
298  reverseLoop = false;
299  break;
300  }
301  else if (index1 == f.rcIndex(index0))
302  {
303  masterIndex = facei;
304  reverseLoop = true;
305  break;
306  }
307  }
308  }
309  }
310 
311  if (masterIndex == -1)
312  {
313  writeOBJ(fp, mesh_.time().path()/"facesToBeMerged.obj");
315  << "Problem" << abort(FatalError);
316  }
317 
318 
319  // Modify the master face.
320  // ~~~~~~~~~~~~~~~~~~~~~~~
321 
322  // Modify first face.
323  label facei = faceLabels[masterIndex];
324 
325  label own = mesh_.faceOwner()[facei];
326 
327  if (cellRegion[own] != -1)
328  {
329  own = cellRegionMaster[cellRegion[own]];
330  }
331 
332  label patchID, zoneID, zoneFlip;
333 
334  getFaceInfo(facei, patchID, zoneID, zoneFlip);
335 
336  label nei = -1;
337 
338  if (mesh_.isInternalFace(facei))
339  {
340  nei = mesh_.faceNeighbour()[facei];
341 
342  if (cellRegion[nei] != -1)
343  {
344  nei = cellRegionMaster[cellRegion[nei]];
345  }
346  }
347 
348 
349  DynamicList<label> faceVerts(edgeLoop.size());
350 
351  forAll(edgeLoop, i)
352  {
353  label pointi = fp.meshPoints()[edgeLoop[i]];
354 
355  if (pointsToRemove.found(pointi))
356  {
357  //Pout<< "**Removing point " << pointi << " from "
358  // << edgeLoop << endl;
359  }
360  else
361  {
362  faceVerts.append(pointi);
363  }
364  }
365 
366  face mergedFace;
367  mergedFace.transfer(faceVerts);
368 
369  if (reverseLoop)
370  {
371  reverse(mergedFace);
372  }
373 
374  //{
375  // Pout<< "Modifying masterface " << facei
376  // << " from faces:" << faceLabels
377  // << " old verts:" << UIndirectList<face>(mesh_.faces(), faceLabels)
378  // << " for new verts:"
379  // << mergedFace
380  // << " possibly new owner " << own
381  // << " or new nei " << nei
382  // << endl;
383  //}
384 
385  modFace
386  (
387  mergedFace, // modified face
388  facei, // label of face being modified
389  own, // owner
390  nei, // neighbour
391  false, // face flip
392  patchID, // patch for face
393  false, // remove from zone
394  zoneID, // zone for face
395  zoneFlip, // face flip in zone
396 
397  meshMod
398  );
399 
400 
401  // Remove all but master face.
402  forAll(faceLabels, patchFacei)
403  {
404  if (patchFacei != masterIndex)
405  {
406  //Pout<< "Removing face " << faceLabels[patchFacei] << endl;
407 
408  meshMod.setAction(polyRemoveFace(faceLabels[patchFacei], facei));
409  }
410  }
411 }
412 
413 
414 // Get patch, zone info for facei
415 void Foam::removeFaces::getFaceInfo
416 (
417  const label facei,
418 
419  label& patchID,
420  label& zoneID,
421  label& zoneFlip
422 ) const
423 {
424  patchID = -1;
425 
426  if (!mesh_.isInternalFace(facei))
427  {
428  patchID = mesh_.boundaryMesh().whichPatch(facei);
429  }
430 
431  zoneID = mesh_.faceZones().whichZone(facei);
432 
433  zoneFlip = false;
434 
435  if (zoneID >= 0)
436  {
437  const faceZone& fZone = mesh_.faceZones()[zoneID];
438 
439  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
440  }
441 }
442 
443 
444 // Return face with all pointsToRemove removed.
445 Foam::face Foam::removeFaces::filterFace
446 (
447  const labelHashSet& pointsToRemove,
448  const label facei
449 ) const
450 {
451  const face& f = mesh_.faces()[facei];
452 
453  labelList newFace(f.size(), -1);
454 
455  label newFp = 0;
456 
457  forAll(f, fp)
458  {
459  label vertI = f[fp];
460 
461  if (!pointsToRemove.found(vertI))
462  {
463  newFace[newFp++] = vertI;
464  }
465  }
466 
467  newFace.setSize(newFp);
468 
469  return face(newFace);
470 }
471 
472 
473 // Wrapper for meshMod.modifyFace. Reverses face if own>nei.
474 void Foam::removeFaces::modFace
475 (
476  const face& f,
477  const label masterFaceID,
478  const label own,
479  const label nei,
480  const bool flipFaceFlux,
481  const label newPatchID,
482  const bool removeFromZone,
483  const label zoneID,
484  const bool zoneFlip,
485 
486  polyTopoChange& meshMod
487 ) const
488 {
489  if ((nei == -1) || (own < nei))
490  {
491 // if (debug)
492 // {
493 // Pout<< "ModifyFace (unreversed) :"
494 // << " facei:" << masterFaceID
495 // << " f:" << f
496 // << " own:" << own
497 // << " nei:" << nei
498 // << " flipFaceFlux:" << flipFaceFlux
499 // << " newPatchID:" << newPatchID
500 // << " removeFromZone:" << removeFromZone
501 // << " zoneID:" << zoneID
502 // << " zoneFlip:" << zoneFlip
503 // << endl;
504 // }
505 
506  meshMod.setAction
507  (
508  polyModifyFace
509  (
510  f, // modified face
511  masterFaceID, // label of face being modified
512  own, // owner
513  nei, // neighbour
514  flipFaceFlux, // face flip
515  newPatchID, // patch for face
516  removeFromZone, // remove from zone
517  zoneID, // zone for face
518  zoneFlip // face flip in zone
519  )
520  );
521  }
522  else
523  {
524 // if (debug)
525 // {
526 // Pout<< "ModifyFace (!reversed) :"
527 // << " facei:" << masterFaceID
528 // << " f:" << f.reverseFace()
529 // << " own:" << nei
530 // << " nei:" << own
531 // << " flipFaceFlux:" << flipFaceFlux
532 // << " newPatchID:" << newPatchID
533 // << " removeFromZone:" << removeFromZone
534 // << " zoneID:" << zoneID
535 // << " zoneFlip:" << zoneFlip
536 // << endl;
537 // }
538 
539  meshMod.setAction
540  (
541  polyModifyFace
542  (
543  f.reverseFace(),// modified face
544  masterFaceID, // label of face being modified
545  nei, // owner
546  own, // neighbour
547  flipFaceFlux, // face flip
548  newPatchID, // patch for face
549  removeFromZone, // remove from zone
550  zoneID, // zone for face
551  zoneFlip // face flip in zone
552  )
553  );
554  }
555 }
556 
557 
558 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
559 
560 // Construct from mesh
561 Foam::removeFaces::removeFaces
562 (
563  const polyMesh& mesh,
564  const scalar minCos
565 )
566 :
567  mesh_(mesh),
568  minCos_(minCos)
569 {}
570 
571 
572 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
573 
574 // Removing face connects cells. This function works out a consistent set of
575 // cell regions.
576 // - returns faces to remove. Can be extended with additional faces
577 // (if owner would become neighbour)
578 // - sets cellRegion to -1 or to region number
579 // - regionMaster contains for every region number a master cell.
581 (
582  const labelList& facesToRemove,
583  labelList& cellRegion,
584  labelList& regionMaster,
585  labelList& newFacesToRemove
586 ) const
587 {
588  const labelList& faceOwner = mesh_.faceOwner();
589  const labelList& faceNeighbour = mesh_.faceNeighbour();
590 
591  cellRegion.setSize(mesh_.nCells());
592  cellRegion = -1;
593 
594  regionMaster.setSize(mesh_.nCells());
595  regionMaster = -1;
596 
597  label nRegions = 0;
598 
599  forAll(facesToRemove, i)
600  {
601  label facei = facesToRemove[i];
602 
603  if (!mesh_.isInternalFace(facei))
604  {
606  << "Not internal face:" << facei << abort(FatalError);
607  }
608 
609 
610  label own = faceOwner[facei];
611  label nei = faceNeighbour[facei];
612 
613  label region0 = cellRegion[own];
614  label region1 = cellRegion[nei];
615 
616  if (region0 == -1)
617  {
618  if (region1 == -1)
619  {
620  // Create new region
621  cellRegion[own] = nRegions;
622  cellRegion[nei] = nRegions;
623 
624  // Make owner (lowest numbered!) the master of the region
625  regionMaster[nRegions] = own;
626  nRegions++;
627  }
628  else
629  {
630  // Add owner to neighbour region
631  cellRegion[own] = region1;
632  // See if owner becomes the master of the region
633  regionMaster[region1] = min(own, regionMaster[region1]);
634  }
635  }
636  else
637  {
638  if (region1 == -1)
639  {
640  // Add neighbour to owner region
641  cellRegion[nei] = region0;
642  // nei is higher numbered than own so guaranteed not lower
643  // than master of region0.
644  }
645  else if (region0 != region1)
646  {
647  // Both have regions. Keep lowest numbered region and master.
648  label freedRegion = -1;
649  label keptRegion = -1;
650 
651  if (region0 < region1)
652  {
653  changeCellRegion
654  (
655  nei,
656  region1, // old region
657  region0, // new region
658  cellRegion
659  );
660 
661  keptRegion = region0;
662  freedRegion = region1;
663  }
664  else if (region1 < region0)
665  {
666  changeCellRegion
667  (
668  own,
669  region0, // old region
670  region1, // new region
671  cellRegion
672  );
673 
674  keptRegion = region1;
675  freedRegion = region0;
676  }
677 
678  label master0 = regionMaster[region0];
679  label master1 = regionMaster[region1];
680 
681  regionMaster[freedRegion] = -1;
682  regionMaster[keptRegion] = min(master0, master1);
683  }
684  }
685  }
686 
687  regionMaster.setSize(nRegions);
688 
689 
690  // Various checks
691  // - master is lowest numbered in any region
692  // - regions have more than 1 cell
693  {
694  labelList nCells(regionMaster.size(), 0);
695 
696  forAll(cellRegion, celli)
697  {
698  label r = cellRegion[celli];
699 
700  if (r != -1)
701  {
702  nCells[r]++;
703 
704  if (celli < regionMaster[r])
705  {
707  << "Not lowest numbered : cell:" << celli
708  << " region:" << r
709  << " regionmaster:" << regionMaster[r]
710  << abort(FatalError);
711  }
712  }
713  }
714 
715  forAll(nCells, region)
716  {
717  if (nCells[region] == 1)
718  {
720  << "Region " << region
721  << " has only " << nCells[region] << " cells in it"
722  << abort(FatalError);
723  }
724  }
725  }
726 
727 
728  // Count number of used regions
729  label nUsedRegions = 0;
730 
731  forAll(regionMaster, i)
732  {
733  if (regionMaster[i] != -1)
734  {
735  nUsedRegions++;
736  }
737  }
738 
739  // Recreate facesToRemove to be consistent with the cellRegions.
740  DynamicList<label> allFacesToRemove(facesToRemove.size());
741 
742  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
743  {
744  label own = faceOwner[facei];
745  label nei = faceNeighbour[facei];
746 
747  if (cellRegion[own] != -1 && cellRegion[own] == cellRegion[nei])
748  {
749  // Both will become the same cell so add face to list of faces
750  // to be removed.
751  allFacesToRemove.append(facei);
752  }
753  }
754 
755  newFacesToRemove.transfer(allFacesToRemove);
756 
757  return nUsedRegions;
758 }
759 
760 
762 (
763  const labelList& faceLabels,
764  const labelList& cellRegion,
765  const labelList& cellRegionMaster,
766  polyTopoChange& meshMod
767 ) const
768 {
769  if (debug)
770  {
771  faceSet facesToRemove(mesh_, "facesToRemove", faceLabels);
772  Pout<< "Writing faces to remove to faceSet " << facesToRemove.name()
773  << endl;
774  facesToRemove.write();
775  }
776 
777  // Make map of all faces to be removed
778  boolList removedFace(mesh_.nFaces(), false);
779 
780  forAll(faceLabels, i)
781  {
782  label facei = faceLabels[i];
783 
784  if (!mesh_.isInternalFace(facei))
785  {
787  << "Face to remove is not internal face:" << facei
788  << abort(FatalError);
789  }
790 
791  removedFace[facei] = true;
792  }
793 
794 
795  // Edges to be removed
796  // ~~~~~~~~~~~~~~~~~~~
797 
798 
799  // Edges to remove
800  labelHashSet edgesToRemove(faceLabels.size());
801 
802  // Per face the region it is in. -1 for removed faces, -2 for regions
803  // consisting of single face only.
804  labelList faceRegion(mesh_.nFaces(), -1);
805 
806  // Number of connected face regions
807  label nRegions = 0;
808 
809  // Storage for on-the-fly addressing
812 
813 
814  {
815  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
816 
817  // Usage of edges by non-removed faces.
818  // See below about initialization.
819  labelList nFacesPerEdge(mesh_.nEdges(), -1);
820 
821  // Count usage of edges by non-removed faces.
822  forAll(faceLabels, i)
823  {
824  label facei = faceLabels[i];
825 
826  const labelList& fEdges = mesh_.faceEdges(facei, fe);
827 
828  forAll(fEdges, i)
829  {
830  label edgeI = fEdges[i];
831 
832  if (nFacesPerEdge[edgeI] == -1)
833  {
834  nFacesPerEdge[edgeI] = mesh_.edgeFaces(edgeI, ef).size()-1;
835  }
836  else
837  {
838  nFacesPerEdge[edgeI]--;
839  }
840  }
841  }
842 
843  // Count usage for edges not on faces-to-be-removed.
844  // Note that this only needs to be done for possibly coupled edges
845  // so we could choose to loop only over boundary faces and use faceEdges
846  // of those.
847 
848  forAll(mesh_.edges(), edgeI)
849  {
850  if (nFacesPerEdge[edgeI] == -1)
851  {
852  // Edge not yet handled in loop above so is not used by any
853  // face to be removed.
854 
855  const labelList& eFaces = mesh_.edgeFaces(edgeI, ef);
856 
857  if (eFaces.size() > 2)
858  {
859  nFacesPerEdge[edgeI] = eFaces.size();
860  }
861  else if (eFaces.size() == 2)
862  {
863  // nFacesPerEdge already -1 so do nothing.
864  }
865  else
866  {
867  const edge& e = mesh_.edges()[edgeI];
868 
870  << "Problem : edge has too few face neighbours:"
871  << eFaces << endl
872  << "edge:" << edgeI
873  << " vertices:" << e
874  << " coords:" << mesh_.points()[e[0]]
875  << mesh_.points()[e[1]]
876  << abort(FatalError);
877  }
878  }
879  }
880 
881 
882 
883  if (debug)
884  {
885  OFstream str(mesh_.time().path()/"edgesWithTwoFaces.obj");
886  Pout<< "Dumping edgesWithTwoFaces to " << str.name() << endl;
887  label vertI = 0;
888 
889  forAll(nFacesPerEdge, edgeI)
890  {
891  if (nFacesPerEdge[edgeI] == 2)
892  {
893  // Edge will get removed.
894  const edge& e = mesh_.edges()[edgeI];
895 
896  meshTools::writeOBJ(str, mesh_.points()[e[0]]);
897  vertI++;
898  meshTools::writeOBJ(str, mesh_.points()[e[1]]);
899  vertI++;
900  str<< "l " << vertI-1 << ' ' << vertI << nl;
901  }
902  }
903  }
904 
905 
906  // Now all unaffected edges will have labelMax, all affected edges the
907  // number of unremoved faces.
908 
909  // Filter for edges inbetween two remaining boundary faces that
910  // make too big an angle.
911  forAll(nFacesPerEdge, edgeI)
912  {
913  if (nFacesPerEdge[edgeI] == 2)
914  {
915  // See if they are two boundary faces
916  label f0 = -1;
917  label f1 = -1;
918 
919  const labelList& eFaces = mesh_.edgeFaces(edgeI, ef);
920 
921  forAll(eFaces, i)
922  {
923  label facei = eFaces[i];
924 
925  if (!removedFace[facei] && !mesh_.isInternalFace(facei))
926  {
927  if (f0 == -1)
928  {
929  f0 = facei;
930  }
931  else
932  {
933  f1 = facei;
934  break;
935  }
936  }
937  }
938 
939  if (f0 != -1 && f1 != -1)
940  {
941  // Edge has two boundary faces remaining.
942  // See if should be merged.
943 
944  label patch0 = patches.whichPatch(f0);
945  label patch1 = patches.whichPatch(f1);
946 
947  if (patch0 != patch1)
948  {
949  // Different patches. Do not merge edge.
951  << "not merging faces " << f0 << " and "
952  << f1 << " across patch boundary edge " << edgeI
953  << endl;
954 
955  // Mark so it gets preserved
956  nFacesPerEdge[edgeI] = 3;
957  }
958  else if (minCos_ < 1 && minCos_ > -1)
959  {
960  const polyPatch& pp0 = patches[patch0];
961  const vectorField& n0 = pp0.faceNormals();
962 
963  if
964  (
965  mag
966  (
967  n0[f0 - pp0.start()]
968  & n0[f1 - pp0.start()]
969  )
970  < minCos_
971  )
972  {
974  << "not merging faces " << f0 << " and "
975  << f1 << " across edge " << edgeI
976  << endl;
977 
978  // Angle between two remaining faces too large.
979  // Mark so it gets preserved
980  nFacesPerEdge[edgeI] = 3;
981  }
982  }
983  }
984  else if (f0 != -1 || f1 != -1)
985  {
986  const edge& e = mesh_.edges()[edgeI];
987 
988  // Only found one boundary face. Problem.
990  << "Problem : edge would have one boundary face"
991  << " and one internal face using it." << endl
992  << "Your remove pattern is probably incorrect." << endl
993  << "edge:" << edgeI
994  << " nFaces:" << nFacesPerEdge[edgeI]
995  << " vertices:" << e
996  << " coords:" << mesh_.points()[e[0]]
997  << mesh_.points()[e[1]]
998  << " face0:" << f0
999  << " face1:" << f1
1000  << abort(FatalError);
1001  }
1002  }
1003  }
1004 
1005 
1006 
1007  // Check locally (before synchronizing) for strangeness
1008  forAll(nFacesPerEdge, edgeI)
1009  {
1010  if (nFacesPerEdge[edgeI] == 1)
1011  {
1012  const edge& e = mesh_.edges()[edgeI];
1013 
1015  << "Problem : edge would get 1 face using it only"
1016  << " edge:" << edgeI
1017  << " nFaces:" << nFacesPerEdge[edgeI]
1018  << " vertices:" << e
1019  << " coords:" << mesh_.points()[e[0]]
1020  << ' ' << mesh_.points()[e[1]]
1021  << abort(FatalError);
1022  }
1023  // Could check here for boundary edge with <=1 faces remaining.
1024  }
1025 
1026 
1027  // Synchronize edge usage. This is to make sure that both sides remove
1028  // (or not remove) an edge on the boundary at the same time.
1029  //
1030  // Coupled edges (edge0, edge1 are opposite each other)
1031  // a. edge not on face to be removed, edge has >= 3 faces
1032  // b. ,, edge has 2 faces
1033  // c. edge has >= 3 remaining faces
1034  // d. edge has 2 remaining faces (assume angle>minCos already handled)
1035  //
1036  // - a + a: do not remove edge
1037  // - a + b: do not remove edge
1038  // - a + c: do not remove edge
1039  // - a + d: do not remove edge
1040  //
1041  // - b + b: do not remove edge
1042  // - b + c: do not remove edge
1043  // - b + d: remove edge
1044  //
1045  // - c + c: do not remove edge
1046  // - c + d: do not remove edge
1047  // - d + d: remove edge
1048  //
1049  //
1050  // So code situation a. with >= 3
1051  // b. with -1
1052  // c. with >=3
1053  // d. with 2
1054  // then do max and check result.
1055  //
1056  // a+a : max(3,3) = 3. do not remove
1057  // a+b : max(3,-1) = 3. do not remove
1058  // a+c : max(3,3) = 3. do not remove
1059  // a+d : max(3,2) = 3. do not remove
1060  //
1061  // b+b : max(-1,-1) = -1. do not remove
1062  // b+c : max(-1,3) = 3. do not remove
1063  // b+d : max(-1,2) = 2. remove
1064  //
1065  // c+c : max(3,3) = 3. do not remove
1066  // c+d : max(3,2) = 3. do not remove
1067  //
1068  // d+d : max(2,2) = 2. remove
1069 
1071  (
1072  mesh_,
1073  nFacesPerEdge,
1074  maxEqOp<label>(),
1075  labelMin // guaranteed to be overridden by maxEqOp
1076  );
1077 
1078  // Convert to labelHashSet
1079  forAll(nFacesPerEdge, edgeI)
1080  {
1081  if (nFacesPerEdge[edgeI] == 0)
1082  {
1083  // 0: edge not used anymore.
1084  edgesToRemove.insert(edgeI);
1085  }
1086  else if (nFacesPerEdge[edgeI] == 1)
1087  {
1088  // 1: illegal. Tested above.
1089  }
1090  else if (nFacesPerEdge[edgeI] == 2)
1091  {
1092  // 2: merge faces.
1093  edgesToRemove.insert(edgeI);
1094  }
1095  }
1096 
1097  if (debug)
1098  {
1099  OFstream str(mesh_.time().path()/"edgesToRemove.obj");
1100  Pout<< "Dumping edgesToRemove to " << str.name() << endl;
1101  label vertI = 0;
1102 
1103  forAllConstIter(labelHashSet, edgesToRemove, iter)
1104  {
1105  // Edge will get removed.
1106  const edge& e = mesh_.edges()[iter.key()];
1107 
1108  meshTools::writeOBJ(str, mesh_.points()[e[0]]);
1109  vertI++;
1110  meshTools::writeOBJ(str, mesh_.points()[e[1]]);
1111  vertI++;
1112  str<< "l " << vertI-1 << ' ' << vertI << nl;
1113  }
1114  }
1115 
1116 
1117  // Walk to fill faceRegion with faces that will be connected across
1118  // edges that will be removed.
1119 
1120  label startFacei = 0;
1121 
1122  while (true)
1123  {
1124  // Find unset region.
1125  for (; startFacei < mesh_.nFaces(); startFacei++)
1126  {
1127  if (faceRegion[startFacei] == -1 && !removedFace[startFacei])
1128  {
1129  break;
1130  }
1131  }
1132 
1133  if (startFacei == mesh_.nFaces())
1134  {
1135  break;
1136  }
1137 
1138  // Start walking face-edge-face, crossing edges that will get
1139  // removed. Every thus connected region will get single region
1140  // number.
1141  label nRegion = changeFaceRegion
1142  (
1143  cellRegion,
1144  removedFace,
1145  nFacesPerEdge,
1146  startFacei,
1147  nRegions,
1148  mesh_.faceEdges(startFacei, fe),
1149  faceRegion
1150  );
1151 
1152  if (nRegion < 1)
1153  {
1154  FatalErrorInFunction << "Problem" << abort(FatalError);
1155  }
1156  else if (nRegion == 1)
1157  {
1158  // Reset face to be single region.
1159  faceRegion[startFacei] = -2;
1160  }
1161  else
1162  {
1163  nRegions++;
1164  }
1165  }
1166 
1167 
1168  // Check we're deciding the same on both sides. Since the regioning
1169  // is done based on nFacesPerEdge (which is synced) this should
1170  // indeed be the case.
1171 
1172  labelList nbrFaceRegion(faceRegion);
1173 
1175  (
1176  mesh_,
1177  nbrFaceRegion
1178  );
1179 
1180  labelList toNbrRegion(nRegions, -1);
1181 
1182  for
1183  (
1184  label facei = mesh_.nInternalFaces();
1185  facei < mesh_.nFaces();
1186  facei++
1187  )
1188  {
1189  // Get the neighbouring region.
1190  label nbrRegion = nbrFaceRegion[facei];
1191  label myRegion = faceRegion[facei];
1192 
1193  if (myRegion <= -1 || nbrRegion <= -1)
1194  {
1195  if (nbrRegion != myRegion)
1196  {
1198  << "Inconsistent face region across coupled patches."
1199  << endl
1200  << "This side has for facei:" << facei
1201  << " region:" << myRegion << endl
1202  << "The other side has region:" << nbrRegion
1203  << endl
1204  << "(region -1 means face is to be deleted)"
1205  << abort(FatalError);
1206  }
1207  }
1208  else if (toNbrRegion[myRegion] == -1)
1209  {
1210  // First visit of region. Store correspondence.
1211  toNbrRegion[myRegion] = nbrRegion;
1212  }
1213  else
1214  {
1215  // Second visit of this region.
1216  if (toNbrRegion[myRegion] != nbrRegion)
1217  {
1219  << "Inconsistent face region across coupled patches."
1220  << endl
1221  << "This side has for facei:" << facei
1222  << " region:" << myRegion
1223  << " with coupled neighbouring regions:"
1224  << toNbrRegion[myRegion] << " and "
1225  << nbrRegion
1226  << abort(FatalError);
1227  }
1228  }
1229  }
1230  }
1231 
1232  //if (debug)
1233  //{
1234  // labelListList regionToFaces(invertOneToMany(nRegions, faceRegion));
1235  //
1236  // forAll(regionToFaces, regionI)
1237  // {
1238  // Pout<< " " << regionI << " faces:" << regionToFaces[regionI]
1239  // << endl;
1240  // }
1241  //}
1242 
1243 
1244  // Points to be removed
1245  // ~~~~~~~~~~~~~~~~~~~~
1246 
1247  labelHashSet pointsToRemove(4*faceLabels.size());
1248 
1249 
1250  // Per point count the number of unremoved edges. Store the ones that
1251  // are only used by 2 unremoved edges.
1252  {
1253  // Usage of points by non-removed edges.
1254  labelList nEdgesPerPoint(mesh_.nPoints());
1255 
1256  const labelListList& pointEdges = mesh_.pointEdges();
1257 
1258  forAll(pointEdges, pointi)
1259  {
1260  nEdgesPerPoint[pointi] = pointEdges[pointi].size();
1261  }
1262 
1263  forAllConstIter(labelHashSet, edgesToRemove, iter)
1264  {
1265  // Edge will get removed.
1266  const edge& e = mesh_.edges()[iter.key()];
1267 
1268  forAll(e, i)
1269  {
1270  nEdgesPerPoint[e[i]]--;
1271  }
1272  }
1273 
1274  // Check locally (before synchronizing) for strangeness
1275  forAll(nEdgesPerPoint, pointi)
1276  {
1277  if (nEdgesPerPoint[pointi] == 1)
1278  {
1280  << "Problem : point would get 1 edge using it only."
1281  << " pointi:" << pointi
1282  << " coord:" << mesh_.points()[pointi]
1283  << abort(FatalError);
1284  }
1285  }
1286 
1287  // Synchronize point usage. This is to make sure that both sides remove
1288  // (or not remove) a point on the boundary at the same time.
1290  (
1291  mesh_,
1292  nEdgesPerPoint,
1293  maxEqOp<label>(),
1294  labelMin
1295  );
1296 
1297  forAll(nEdgesPerPoint, pointi)
1298  {
1299  if (nEdgesPerPoint[pointi] == 0)
1300  {
1301  pointsToRemove.insert(pointi);
1302  }
1303  else if (nEdgesPerPoint[pointi] == 1)
1304  {
1305  // Already checked before
1306  }
1307  else if (nEdgesPerPoint[pointi] == 2)
1308  {
1309  // Remove point and merge edges.
1310  pointsToRemove.insert(pointi);
1311  }
1312  }
1313  }
1314 
1315 
1316  if (debug)
1317  {
1318  OFstream str(mesh_.time().path()/"pointsToRemove.obj");
1319  Pout<< "Dumping pointsToRemove to " << str.name() << endl;
1320 
1321  forAllConstIter(labelHashSet, pointsToRemove, iter)
1322  {
1323  meshTools::writeOBJ(str, mesh_.points()[iter.key()]);
1324  }
1325  }
1326 
1327 
1328  // All faces affected in any way
1329  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1330 
1331  // Get all faces affected in any way by removal of points/edges/faces/cells
1332  boolList affectedFace
1333  (
1334  getFacesAffected
1335  (
1336  cellRegion,
1337  cellRegionMaster,
1338  faceLabels,
1339  edgesToRemove,
1340  pointsToRemove
1341  )
1342  );
1343 
1344  //
1345  // Now we know
1346  // - faceLabels : faces to remove (sync since no boundary faces)
1347  // - cellRegion/Master : cells to remove (sync since cells)
1348  // - pointsToRemove : points to remove (sync)
1349  // - faceRegion : connected face region of faces to be merged (sync)
1350  // - affectedFace : faces with points removed and/or owner/neighbour
1351  // changed (non sync)
1352 
1353 
1354  // Start modifying mesh and keep track of faces changed.
1355 
1356 
1357  // Do all removals
1358  // ~~~~~~~~~~~~~~~
1359 
1360  // Remove split faces.
1361  forAll(faceLabels, labelI)
1362  {
1363  label facei = faceLabels[labelI];
1364 
1365  // Remove face if not yet uptodate (which is never; but want to be
1366  // consistent with rest of face removals/modifications)
1367  if (affectedFace[facei])
1368  {
1369  affectedFace[facei] = false;
1370 
1371  meshMod.setAction(polyRemoveFace(facei, -1));
1372  }
1373  }
1374 
1375 
1376  // Remove points.
1377  forAllConstIter(labelHashSet, pointsToRemove, iter)
1378  {
1379  label pointi = iter.key();
1380 
1381  meshMod.setAction(polyRemovePoint(pointi, -1));
1382  }
1383 
1384 
1385  // Remove cells.
1386  forAll(cellRegion, celli)
1387  {
1388  label region = cellRegion[celli];
1389 
1390  if (region != -1 && (celli != cellRegionMaster[region]))
1391  {
1392  meshMod.setAction(polyRemoveCell(celli, cellRegionMaster[region]));
1393  }
1394  }
1395 
1396 
1397 
1398  // Merge faces across edges to be merged
1399  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1400 
1401  // Invert faceRegion so we get region to faces.
1402  {
1403  labelListList regionToFaces(invertOneToMany(nRegions, faceRegion));
1404 
1405  forAll(regionToFaces, regionI)
1406  {
1407  const labelList& rFaces = regionToFaces[regionI];
1408 
1409  if (rFaces.size() <= 1)
1410  {
1412  << "Region:" << regionI
1413  << " contains only faces " << rFaces
1414  << abort(FatalError);
1415  }
1416 
1417  // rFaces[0] is master, rest gets removed.
1418  mergeFaces
1419  (
1420  cellRegion,
1421  cellRegionMaster,
1422  pointsToRemove,
1423  rFaces,
1424  meshMod
1425  );
1426 
1427  forAll(rFaces, i)
1428  {
1429  affectedFace[rFaces[i]] = false;
1430  }
1431  }
1432  }
1433 
1434 
1435  // Remaining affected faces
1436  // ~~~~~~~~~~~~~~~~~~~~~~~~
1437 
1438  // Check if any remaining faces have not been updated for new slave/master
1439  // or points removed.
1440  forAll(affectedFace, facei)
1441  {
1442  if (affectedFace[facei])
1443  {
1444  affectedFace[facei] = false;
1445 
1446  face f(filterFace(pointsToRemove, facei));
1447 
1448  label own = mesh_.faceOwner()[facei];
1449 
1450  if (cellRegion[own] != -1)
1451  {
1452  own = cellRegionMaster[cellRegion[own]];
1453  }
1454 
1455  label patchID, zoneID, zoneFlip;
1456 
1457  getFaceInfo(facei, patchID, zoneID, zoneFlip);
1458 
1459  label nei = -1;
1460 
1461  if (mesh_.isInternalFace(facei))
1462  {
1463  nei = mesh_.faceNeighbour()[facei];
1464 
1465  if (cellRegion[nei] != -1)
1466  {
1467  nei = cellRegionMaster[cellRegion[nei]];
1468  }
1469  }
1470 
1471 // if (debug)
1472 // {
1473 // Pout<< "Modifying " << facei
1474 // << " old verts:" << mesh_.faces()[facei]
1475 // << " for new verts:" << f
1476 // << " or for new owner " << own << " or for new nei "
1477 // << nei
1478 // << endl;
1479 // }
1480 
1481  modFace
1482  (
1483  f, // modified face
1484  facei, // label of face being modified
1485  own, // owner
1486  nei, // neighbour
1487  false, // face flip
1488  patchID, // patch for face
1489  false, // remove from zone
1490  zoneID, // zone for face
1491  zoneFlip, // face flip in zone
1492 
1493  meshMod
1494  );
1495  }
1496  }
1497 }
1498 
1499 
1500 // ************************************************************************* //
Class containing data for face removal.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A list of face labels.
Definition: faceSet.H:48
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Output to file stream.
Definition: OFstream.H:81
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
patches[0]
PrimitivePatch< face, IndirectList, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
scalar f1
Definition: createFields.H:28
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Class containing data for cell removal.
static const labelSphericalTensor labelI(1)
Identity labelTensor.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
label compatibleRemoves(const labelList &inPiercedFaces, labelList &cellRegion, labelList &cellRegionMaster, labelList &outPiercedFaces) const
Find faces including those with cells which have the same mastercell.
Definition: removeFaces.C:581
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh edges.
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:97
List< label > labelList
A List of labels.
Definition: labelList.H:56
static void swapFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled face values.
Definition: syncTools.H:463
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::polyBoundaryMesh.
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
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:295
const Field< PointType > & faceNormals() const
Return face normals for patch.
#define WarningInFunction
Report a warning using Foam::Warning.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Direct mesh changes based on v1.3 polyTopoChange syntax.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
virtual bool write() const
Write using setting from DB.
void setRefinement(const labelList &piercedFaces, const labelList &cellRegion, const labelList &cellRegionMaster, polyTopoChange &) const
Play commands into polyTopoChange to remove faces.
Definition: removeFaces.C:762
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:365
Class containing data for point removal.
static const label labelMin
Definition: label.H:61
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.