edgeCollapser.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 "edgeCollapser.H"
27 #include "polyMesh.H"
28 #include "polyTopoChange.H"
29 #include "globalMeshData.H"
30 #include "syncTools.H"
31 #include "PointEdgeWave.H"
32 #include "globalIndex.H"
33 #include "removePoints.H"
34 #include "meshCheck.H"
35 #include "OFstream.H"
36 
37 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
42 }
43 
44 
46 (
47  const polyMesh& mesh,
48  const dictionary& meshQualityDict
49 )
50 {
51  labelHashSet badFaces(mesh.nFaces()/100);
52  DynamicList<label> checkFaces(mesh.nFaces());
53 
54  const vectorField& fAreas = mesh.faceAreas();
55 
56  scalar faceAreaLimit = small;
57 
58  forAll(fAreas, fI)
59  {
60  if (mag(fAreas[fI]) > faceAreaLimit)
61  {
62  checkFaces.append(fI);
63  }
64  }
65 
66  Info<< endl;
67 
69  (
70  false,
71  mesh,
72  meshQualityDict,
73  checkFaces,
74  badFaces
75  );
76 
77  return badFaces;
78 }
79 
80 
82 (
83  const polyMesh& mesh,
84  const dictionary& meshQualityDict,
85  PackedBoolList& isErrorPoint
86 )
87 {
89  (
90  mesh,
91  meshQualityDict
92  );
93 
94  label nBadFaces = returnReduce(badFaces.size(), sumOp<label>());
95 
96  forAllConstIter(labelHashSet, badFaces, iter)
97  {
98  const face& f = mesh.faces()[iter.key()];
99 
100  forAll(f, pI)
101  {
102  isErrorPoint[f[pI]] = true;
103  }
104  }
105 
107  (
108  mesh,
109  isErrorPoint,
111  0
112  );
113 
114  return nBadFaces;
115 }
116 
117 
118 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
119 
120 Foam::labelList Foam::edgeCollapser::edgesFromPoints
121 (
122  const label& facei,
123  const labelList& pointLabels
124 ) const
125 {
126  labelList edgeLabels(pointLabels.size() - 1, -1);
127 
128  const labelList& faceEdges = mesh_.faceEdges()[facei];
129  const edgeList& edges = mesh_.edges();
130 
131  label count = 0;
132 
133  forAll(faceEdges, eI)
134  {
135  const label edgeI = faceEdges[eI];
136  const edge& e = edges[edgeI];
137 
138  label match = 0;
139 
140  forAll(pointLabels, pI)
141  {
142  if (e[0] == pointLabels[pI])
143  {
144  match++;
145  }
146 
147  if (e[1] == pointLabels[pI])
148  {
149  match++;
150  }
151  }
152 
153  if (match == 2)
154  {
155  // Edge found
156  edgeLabels[count++] = edgeI;
157  }
158  }
159 
160  if (count != edgeLabels.size())
161  {
162  edgeLabels.setSize(count);
163  }
164 
165  return edgeLabels;
166 }
167 
168 
169 void Foam::edgeCollapser::collapseToEdge
170 (
171  const label facei,
172  const pointField& pts,
173  const labelList& pointPriority,
174  const vector& collapseAxis,
175  const point& fC,
176  const labelList& facePtsNeg,
177  const labelList& facePtsPos,
178  const scalarList& dNeg,
179  const scalarList& dPos,
180  const scalar dShift,
181  PackedBoolList& collapseEdge,
182  Map<point>& collapsePointToLocation
183 ) const
184 {
185  // Negative half
186 
187  Foam::point collapseToPtA(great, great, great);
188  // collapseAxis*(sum(dNeg)/dNeg.size() - dShift) + fC;
189 
190  label maxPriority = labelMin;
191  DynamicList<label> maxPriorityPts(max(dNeg.size(), dPos.size()));
192 
193  forAll(facePtsNeg, fPtI)
194  {
195  const label facePointi = facePtsNeg[fPtI];
196  const label facePtPriority = pointPriority[facePointi];
197 
198  if (facePtPriority > maxPriority)
199  {
200  maxPriority = facePtPriority;
201  maxPriorityPts.clear();
202  maxPriorityPts.append(facePointi);
203  }
204  else if (facePtPriority == maxPriority)
205  {
206  maxPriorityPts.append(facePointi);
207  }
208  }
209 
210  if (!maxPriorityPts.empty())
211  {
212  Foam::point averagePt(Zero);
213 
214  forAll(maxPriorityPts, ptI)
215  {
216  averagePt += pts[maxPriorityPts[ptI]];
217  }
218 
219  collapseToPtA = averagePt/maxPriorityPts.size();
220 // collapseToPtA = pts[maxPriorityPts.first()];
221  }
222 
223  maxPriority = labelMin;
224  maxPriorityPts.clear();
225 
226  labelList faceEdgesNeg = edgesFromPoints(facei, facePtsNeg);
227 
228  forAll(faceEdgesNeg, edgeI)
229  {
230  collapseEdge[faceEdgesNeg[edgeI]] = true;
231  }
232 
233  forAll(facePtsNeg, pI)
234  {
235  collapsePointToLocation.set(facePtsNeg[pI], collapseToPtA);
236  }
237 
238 
239  // Positive half
240  Foam::point collapseToPtB(great, great, great);
241 // = collapseAxis*(sum(dPos)/dPos.size() - dShift) + fC;
242 
243  forAll(facePtsPos, fPtI)
244  {
245  const label facePointi = facePtsPos[fPtI];
246  const label facePtPriority = pointPriority[facePointi];
247 
248  if (facePtPriority > maxPriority)
249  {
250  maxPriority = facePtPriority;
251  maxPriorityPts.clear();
252  maxPriorityPts.append(facePointi);
253  }
254  else if (facePtPriority == maxPriority)
255  {
256  maxPriorityPts.append(facePointi);
257  }
258  }
259 
260  if (!maxPriorityPts.empty())
261  {
262  Foam::point averagePt(Zero);
263 
264  forAll(maxPriorityPts, ptI)
265  {
266  averagePt += pts[maxPriorityPts[ptI]];
267  }
268 
269  collapseToPtB = averagePt/maxPriorityPts.size();
270 // collapseToPtB = pts[maxPriorityPts.first()];
271  }
272 
273  labelList faceEdgesPos = edgesFromPoints(facei, facePtsPos);
274 
275  forAll(faceEdgesPos, edgeI)
276  {
277  collapseEdge[faceEdgesPos[edgeI]] = true;
278  }
279 
280  forAll(facePtsPos, pI)
281  {
282  collapsePointToLocation.set(facePtsPos[pI], collapseToPtB);
283  }
284 }
285 
286 
287 void Foam::edgeCollapser::collapseToPoint
288 (
289  const label& facei,
290  const pointField& pts,
291  const labelList& pointPriority,
292  const point& fC,
293  const labelList& facePts,
294  PackedBoolList& collapseEdge,
295  Map<point>& collapsePointToLocation
296 ) const
297 {
298  const face& f = mesh_.faces()[facei];
299 
300  Foam::point collapseToPt = fC;
301 
302  label maxPriority = labelMin;
303  DynamicList<label> maxPriorityPts(f.size());
304 
305  forAll(facePts, fPtI)
306  {
307  const label facePointi = facePts[fPtI];
308  const label facePtPriority = pointPriority[facePointi];
309 
310  if (facePtPriority > maxPriority)
311  {
312  maxPriority = facePtPriority;
313  maxPriorityPts.clear();
314  maxPriorityPts.append(facePointi);
315  }
316  else if (facePtPriority == maxPriority)
317  {
318  maxPriorityPts.append(facePointi);
319  }
320  }
321 
322  if (!maxPriorityPts.empty())
323  {
324  Foam::point averagePt(Zero);
325 
326  forAll(maxPriorityPts, ptI)
327  {
328  averagePt += pts[maxPriorityPts[ptI]];
329  }
330 
331  collapseToPt = averagePt/maxPriorityPts.size();
332 
333 // collapseToPt = pts[maxPriorityPts.first()];
334  }
335 
336 // DynamicList<label> faceBoundaryPts(f.size());
337 // DynamicList<label> faceFeaturePts(f.size());
338 //
339 // forAll(facePts, fPtI)
340 // {
341 // if (pointPriority[facePts[fPtI]] == 1)
342 // {
343 // faceFeaturePts.append(facePts[fPtI]);
344 // }
345 // else if (pointPriority[facePts[fPtI]] == 0)
346 // {
347 // faceBoundaryPts.append(facePts[fPtI]);
348 // }
349 // }
350 //
351 // if (!faceBoundaryPts.empty() || !faceFeaturePts.empty())
352 // {
353 // if (!faceFeaturePts.empty())
354 // {
355 // collapseToPt = pts[faceFeaturePts.first()];
356 // }
357 // else if (faceBoundaryPts.size() == 2)
358 // {
359 // collapseToPt =
360 // 0.5
361 // *(
362 // pts[faceBoundaryPts[0]]
363 // + pts[faceBoundaryPts[1]]
364 // );
365 // }
366 // else if (faceBoundaryPts.size() <= f.size())
367 // {
368 // face bFace(faceBoundaryPts);
369 //
370 // collapseToPt = bFace.centre(pts);
371 // }
372 // }
373 
374  const labelList& faceEdges = mesh_.faceEdges()[facei];
375 
376  forAll(faceEdges, eI)
377  {
378  const label edgeI = faceEdges[eI];
379  collapseEdge[edgeI] = true;
380  }
381 
382  forAll(f, pI)
383  {
384  collapsePointToLocation.set(f[pI], collapseToPt);
385  }
386 }
387 
388 
389 void Foam::edgeCollapser::faceCollapseAxisAndAspectRatio
390 (
391  const face& f,
392  const point& fC,
393  vector& collapseAxis,
394  scalar& aspectRatio
395 ) const
396 {
397  const pointField& pts = mesh_.points();
398 
399  tensor J = f.inertia(pts, fC);
400 
401  // Find the dominant collapse direction by finding the eigenvector
402  // that corresponds to the normal direction, discarding it. The
403  // eigenvector corresponding to the smaller of the two remaining
404  // eigenvalues is the dominant axis in a high aspect ratio face.
405 
406  scalar magJ = mag(J);
407 
408  scalar detJ = small;
409 
410  if (magJ > vSmall)
411  {
412  // Normalise inertia tensor to remove problems with small values
413 
414  J /= mag(J);
415  // J /= cmptMax(J);
416  // J /= max(eigenValues(J).x(), small);
417 
418  // Calculating determinant, including stabilisation for zero or
419  // small negative values
420 
421  detJ = max(det(J), small);
422  }
423 
424  if (detJ < 1e-5)
425  {
426  collapseAxis = f.edges()[longestEdge(f, pts)].vec(pts);
427 
428  // It is possible that all the points of a face are the same
429  if (magSqr(collapseAxis) > vSmall)
430  {
431  collapseAxis /= mag(collapseAxis);
432  }
433 
434  // Empirical correlation for high aspect ratio faces
435 
436  aspectRatio = Foam::sqrt(0.35/detJ);
437  }
438  else
439  {
440  vector eVals = eigenValues(J);
441 
442  if (mag(eVals.y() - eVals.x()) < 100*small)
443  {
444  // First two eigenvalues are the same: i.e. a square face
445 
446  // Cannot necessarily determine linearly independent
447  // eigenvectors, or any at all, use longest edge direction.
448 
449  collapseAxis = f.edges()[longestEdge(f, pts)].vec(pts);
450 
451  collapseAxis /= mag(collapseAxis);
452 
453  aspectRatio = 1.0;
454  }
455  else
456  {
457  // The maximum eigenvalue (z()) must be the direction of the
458  // normal, as it has the greatest value. The minimum eigenvalue
459  // is the dominant collapse axis for high aspect ratio faces.
460 
461  collapseAxis = eigenVectors(J, eVals).x();
462 
463  // The inertia calculation describes the mass distribution as a
464  // function of distance squared to the axis, so the square root of
465  // the ratio of face-plane moments gives a good indication of the
466  // aspect ratio.
467 
468  aspectRatio = Foam::sqrt(eVals.y()/max(eVals.x(), small));
469  }
470  }
471 }
472 
473 
474 Foam::scalarField Foam::edgeCollapser::calcTargetFaceSizes() const
475 {
476  scalarField targetFaceSizes(mesh_.nFaces(), -1);
477 
478  const scalarField& V = mesh_.cellVolumes();
479  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
480 
481  const labelList& cellOwner = mesh_.faceOwner();
482  const labelList& cellNeighbour = mesh_.faceNeighbour();
483 
484  const label nBoundaryFaces = mesh_.nFaces() - mesh_.nInternalFaces();
485 
486  // Calculate face size from cell volumes for internal faces
487  for (label intFacei = 0; intFacei < mesh_.nInternalFaces(); ++intFacei)
488  {
489  const scalar cellOwnerVol = max(0.0, V[cellOwner[intFacei]]);
490  const scalar cellNeighbourVol = max(0.0, V[cellNeighbour[intFacei]]);
491 
492  scalar targetFaceSizeA = Foam::pow(cellOwnerVol, 1.0/3.0);
493  scalar targetFaceSizeB = Foam::pow(cellNeighbourVol, 1.0/3.0);
494 
495  targetFaceSizes[intFacei] = 0.5*(targetFaceSizeA + targetFaceSizeB);
496  }
497 
498  scalarField neiCellVolumes(nBoundaryFaces, -1);
499 
500  // Now do boundary faces
502  {
503  const polyPatch& patch = patches[patchi];
504 
505  label bFacei = patch.start() - mesh_.nInternalFaces();
506 
507  if (patch.coupled())
508  {
509  // Processor boundary face: Need to get the cell volume on the other
510  // processor
511  const labelUList& faceCells = patch.faceCells();
512 
513  forAll(faceCells, facei)
514  {
515  neiCellVolumes[bFacei++] = max(0.0, V[faceCells[facei]]);
516  }
517  }
518  else
519  {
520  // Normal boundary face: Just use owner cell volume to calculate
521  // the target face size
522  forAll(patch, patchFacei)
523  {
524  const label extFacei = patchFacei + patch.start();
525  const scalar cellOwnerVol = max(0.0, V[cellOwner[extFacei]]);
526 
527  targetFaceSizes[extFacei] = Foam::pow(cellOwnerVol, 1.0/3.0);
528  }
529  }
530  }
531 
532  syncTools::swapBoundaryFaceList(mesh_, neiCellVolumes);
533 
535  {
536  const polyPatch& patch = patches[patchi];
537 
538  label bFacei = patch.start() - mesh_.nInternalFaces();
539 
540  if (patch.coupled())
541  {
542  forAll(patch, patchFacei)
543  {
544  const label localFacei = patchFacei + patch.start();
545  const scalar cellOwnerVol = max(0.0, V[cellOwner[localFacei]]);
546  const scalar cellNeighbourVol = neiCellVolumes[bFacei++];
547 
548  scalar targetFaceSizeA = Foam::pow(cellOwnerVol, 1.0/3.0);
549  scalar targetFaceSizeB = Foam::pow(cellNeighbourVol, 1.0/3.0);
550 
551  targetFaceSizes[localFacei]
552  = 0.5*(targetFaceSizeA + targetFaceSizeB);
553  }
554  }
555  }
556 
557  // Returns a characteristic length, not an area
558  return targetFaceSizes;
559 }
560 
561 
562 Foam::edgeCollapser::collapseType Foam::edgeCollapser::collapseFace
563 (
564  const labelList& pointPriority,
565  const face& f,
566  const label facei,
567  const scalar targetFaceSize,
568  PackedBoolList& collapseEdge,
569  Map<point>& collapsePointToLocation,
570  const scalarField& faceFilterFactor
571 ) const
572 {
573  const scalar collapseSizeLimitCoeff = faceFilterFactor[facei];
574 
575  const pointField& pts = mesh_.points();
576 
577  labelList facePts(f);
578 
579  const Foam::point fC = f.centre(pts);
580 
581  const scalar fA = f.mag(pts);
582 
583  vector collapseAxis = Zero;
584  scalar aspectRatio = 1.0;
585 
586  faceCollapseAxisAndAspectRatio(f, fC, collapseAxis, aspectRatio);
587 
588  // The signed distance along the collapse axis passing through the
589  // face centre that each vertex projects to.
590 
591  scalarField d(f.size());
592 
593  forAll(f, fPtI)
594  {
595  const Foam::point& pt = pts[f[fPtI]];
596 
597  d[fPtI] = (collapseAxis & (pt - fC));
598  }
599 
600  // Sort the projected distances and the corresponding vertex
601  // indices along the collapse axis
602 
603  labelList oldToNew;
604 
605  sortedOrder(d, oldToNew);
606 
607  oldToNew = invert(oldToNew.size(), oldToNew);
608 
609  inplaceReorder(oldToNew, d);
610 
611  inplaceReorder(oldToNew, facePts);
612 
613  // Shift the points so that they are relative to the centre of the
614  // collapse line.
615 
616  scalar dShift = -0.5*(d.first() + d.last());
617 
618  d += dShift;
619 
620  // Form two lists, one for each half of the set of points
621  // projected along the collapse axis.
622 
623  // Middle value, index of first entry in the second half
624  label middle = -1;
625 
626  forAll(d, dI)
627  {
628  if (d[dI] > 0)
629  {
630  middle = dI;
631 
632  break;
633  }
634  }
635 
636  if (middle == -1)
637  {
638 // SeriousErrorInFunction
639 // << "middle == -1, " << f << " " << d
640 // << endl;//abort(FatalError);
641 
642  return noCollapse;
643  }
644 
645  // Negative half
646  SubList<scalar> dNeg(d, middle, 0);
647  SubList<label> facePtsNeg(facePts, middle, 0);
648 
649  // Positive half
650  SubList<scalar> dPos(d, d.size() - middle, middle);
651  SubList<label> facePtsPos(facePts, d.size() - middle, middle);
652 
653  // Defining how close to the midpoint (M) of the projected
654  // vertices line a projected vertex (X) can be before making this
655  // an invalid edge collapse
656  //
657  // X---X-g----------------M----X-----------g----X--X
658  //
659  // Only allow a collapse if all projected vertices are outwith
660  // guardFraction (g) of the distance form the face centre to the
661  // furthest vertex in the considered direction
662 
663  if (dNeg.size() == 0 || dPos.size() == 0)
664  {
666  << "All points on one side of face centre, not collapsing."
667  << endl;
668  }
669 
670 // Info<< "Face : " << f << nl
671 // << " Collapse Axis: " << collapseAxis << nl
672 // << " Aspect Ratio : " << aspectRatio << endl;
673 
674  collapseType typeOfCollapse = noCollapse;
675 
676  if (magSqr(collapseAxis) < vSmall)
677  {
678  typeOfCollapse = toPoint;
679  }
680  else if (fA < aspectRatio*sqr(targetFaceSize*collapseSizeLimitCoeff))
681  {
682  if
683  (
684  allowEarlyCollapseToPoint_
685  && (d.last() - d.first())
686  < targetFaceSize
687  *allowEarlyCollapseCoeff_*maxCollapseFaceToPointSideLengthCoeff_
688  )
689  {
690  typeOfCollapse = toPoint;
691  }
692  else if
693  (
694  (dNeg.last() < guardFraction_*dNeg.first())
695  && (dPos.first() > guardFraction_*dPos.last())
696  )
697  {
698  typeOfCollapse = toEdge;
699  }
700  else if
701  (
702  (d.last() - d.first())
703  < targetFaceSize
704  *maxCollapseFaceToPointSideLengthCoeff_
705  )
706  {
707  // If the face can't be collapsed to an edge, and it has a
708  // small enough span, collapse it to a point.
709  typeOfCollapse = toPoint;
710  }
711  }
712 
713  if (typeOfCollapse == toPoint)
714  {
715  collapseToPoint
716  (
717  facei,
718  pts,
719  pointPriority,
720  fC,
721  facePts,
722  collapseEdge,
723  collapsePointToLocation
724  );
725  }
726  else if (typeOfCollapse == toEdge)
727  {
728  collapseToEdge
729  (
730  facei,
731  pts,
732  pointPriority,
733  collapseAxis,
734  fC,
735  facePtsNeg,
736  facePtsPos,
737  dNeg,
738  dPos,
739  dShift,
740  collapseEdge,
741  collapsePointToLocation
742  );
743  }
744 
745  return typeOfCollapse;
746 }
747 
748 
749 Foam::label Foam::edgeCollapser::edgeMaster
750 (
751  const labelList& pointPriority,
752  const edge& e
753 ) const
754 {
755  label masterPoint = -1;
756 
757  const label e0 = e.start();
758  const label e1 = e.end();
759 
760  const label e0Priority = pointPriority[e0];
761  const label e1Priority = pointPriority[e1];
762 
763  if (e0Priority > e1Priority)
764  {
765  masterPoint = e0;
766  }
767  else if (e0Priority < e1Priority)
768  {
769  masterPoint = e1;
770  }
771  else if (e0Priority == e1Priority)
772  {
773  masterPoint = e0;
774  }
775 
776 // // Collapse edge to point with higher priority.
777 // if (pointPriority[e0] >= 0)
778 // {
779 // if (pointPriority[e1] >= 0)
780 // {
781 // // Both points have high priority. Choose one to collapse to.
782 // // Note: should look at feature edges/points!
783 // masterPoint = e0;
784 // }
785 // else
786 // {
787 // masterPoint = e0;
788 // }
789 // }
790 // else
791 // {
792 // if (pointPriority[e1] >= 0)
793 // {
794 // masterPoint = e1;
795 // }
796 // else
797 // {
798 // // None on boundary. Neither is a master.
799 // return -1;
800 // }
801 // }
802 
803  return masterPoint;
804 }
805 
806 
807 void Foam::edgeCollapser::checkBoundaryPointMergeEdges
808 (
809  const label pointi,
810  const label otherPointi,
811  const labelList& pointPriority,
812  Map<point>& collapsePointToLocation
813 ) const
814 {
815  const pointField& points = mesh_.points();
816 
817  const label e0Priority = pointPriority[pointi];
818  const label e1Priority = pointPriority[otherPointi];
819 
820  if (e0Priority > e1Priority)
821  {
822  collapsePointToLocation.set
823  (
824  otherPointi,
825  points[pointi]
826  );
827  }
828  else if (e0Priority < e1Priority)
829  {
830  collapsePointToLocation.set
831  (
832  pointi,
833  points[otherPointi]
834  );
835  }
836  else // e0Priority == e1Priority
837  {
838  collapsePointToLocation.set
839  (
840  pointi,
841  points[otherPointi]
842  );
843 
844 // Foam::point averagePt
845 // (
846 // 0.5*(points[otherPointi] + points[pointi])
847 // );
848 //
849 // collapsePointToLocation.set(pointi, averagePt);
850 // collapsePointToLocation.set(otherPointi, averagePt);
851  }
852 }
853 
854 
855 Foam::label Foam::edgeCollapser::breakStringsAtEdges
856 (
857  PackedBoolList& collapseEdge,
858  List<pointEdgeCollapse>& allPointInfo
859 ) const
860 {
861  const edgeList& edges = mesh_.edges();
862  const labelListList& pointEdges = mesh_.pointEdges();
863 
864  label nUncollapsed = 0;
865 
866  forAll(edges, eI)
867  {
868  const edge& e = edges[eI];
869 
870  const label startCollapseIndex
871  = allPointInfo[e.start()].collapseIndex();
872 
873  if (startCollapseIndex != -1 && startCollapseIndex != -2)
874  {
875  const label endCollapseIndex
876  = allPointInfo[e.end()].collapseIndex();
877 
878  if (!collapseEdge[eI] && startCollapseIndex == endCollapseIndex)
879  {
880  const labelList& ptEdgesStart = pointEdges[e.start()];
881 
882  forAll(ptEdgesStart, ptEdgeI)
883  {
884  const label edgeI = ptEdgesStart[ptEdgeI];
885 
886  const label nbrPointi = edges[edgeI].otherVertex(e.start());
887  const label nbrIndex =
888  allPointInfo[nbrPointi].collapseIndex();
889 
890  if (collapseEdge[edgeI] && nbrIndex == startCollapseIndex)
891  {
892  collapseEdge[edgeI] = false;
893  nUncollapsed++;
894  }
895  }
896  }
897  }
898  }
899 
900  return nUncollapsed;
901 }
902 
903 
904 void Foam::edgeCollapser::determineDuplicatePointsOnFace
905 (
906  const face& f,
907  PackedBoolList& markedPoints,
908  labelHashSet& uniqueCollapses,
909  labelHashSet& duplicateCollapses,
910  List<pointEdgeCollapse>& allPointInfo
911 ) const
912 {
913  uniqueCollapses.clear();
914  duplicateCollapses.clear();
915 
916  forAll(f, fpI)
917  {
918  label index = allPointInfo[f[fpI]].collapseIndex();
919 
920  // Check for consecutive duplicate
921  if (index != allPointInfo[f.prevLabel(fpI)].collapseIndex())
922  {
923  if (!uniqueCollapses.insert(index))
924  {
925  // Failed inserting so duplicate
926  duplicateCollapses.insert(index);
927  }
928  }
929  }
930 
931  // Now duplicateCollapses contains duplicate collapse indices.
932  // Convert to points.
933  forAll(f, fpI)
934  {
935  label index = allPointInfo[f[fpI]].collapseIndex();
936  if (duplicateCollapses.found(index))
937  {
938  markedPoints[f[fpI]] = true;
939  }
940  }
941 }
942 
943 
944 Foam::label Foam::edgeCollapser::countEdgesOnFace
945 (
946  const face& f,
947  List<pointEdgeCollapse>& allPointInfo
948 ) const
949 {
950  label nEdges = 0;
951 
952  forAll(f, fpI)
953  {
954  const label pointi = f[fpI];
955  const label newPointi = allPointInfo[pointi].collapseIndex();
956 
957  if (newPointi == -2)
958  {
959  nEdges++;
960  }
961  else
962  {
963  const label prevPointi = f[f.fcIndex(fpI)];
964  const label prevNewPointi
965  = allPointInfo[prevPointi].collapseIndex();
966 
967  if (newPointi != prevNewPointi)
968  {
969  nEdges++;
970  }
971  }
972  }
973 
974  return nEdges;
975 }
976 
977 
978 bool Foam::edgeCollapser::isFaceCollapsed
979 (
980  const face& f,
981  List<pointEdgeCollapse>& allPointInfo
982 ) const
983 {
984  label nEdges = countEdgesOnFace(f, allPointInfo);
985 
986  // Polygons must have 3 or more edges to be valid
987  if (nEdges < 3)
988  {
989  return true;
990  }
991 
992  return false;
993 }
994 
995 
996 // Create consistent set of collapses.
997 // collapseEdge : per edge:
998 // -1 : do not collapse
999 // 0 : collapse to start
1000 // 1 : collapse to end
1001 // Note: collapseEdge has to be parallel consistent (in orientation)
1002 Foam::label Foam::edgeCollapser::syncCollapse
1003 (
1004  const globalIndex& globalPoints,
1005  const labelList& pointPriority,
1006  const PackedBoolList& collapseEdge,
1007  const Map<point>& collapsePointToLocation,
1008  List<pointEdgeCollapse>& allPointInfo
1009 ) const
1010 {
1011  const edgeList& edges = mesh_.edges();
1012 
1013  label nCollapsed = 0;
1014 
1015  DynamicList<label> initPoints(mesh_.nPoints());
1016  DynamicList<pointEdgeCollapse> initPointInfo(mesh_.nPoints());
1017 
1018  allPointInfo.clear();
1019  allPointInfo.setSize(mesh_.nPoints());
1020 
1021  // Initialise edges to no collapse
1022  List<pointEdgeCollapse> allEdgeInfo
1023  (
1024  mesh_.nEdges(),
1025  pointEdgeCollapse(Zero, -1, -1)
1026  );
1027 
1028  // Mark selected edges for collapse
1029  forAll(edges, edgeI)
1030  {
1031  if (collapseEdge[edgeI])
1032  {
1033  const edge& e = edges[edgeI];
1034 
1035  label masterPointi = e.start();
1036 
1037  // Choose the point on the edge with the highest priority.
1038  if (pointPriority[e.end()] > pointPriority[e.start()])
1039  {
1040  masterPointi = e.end();
1041  }
1042 
1043  label masterPointPriority = pointPriority[masterPointi];
1044 
1045  label index = globalPoints.toGlobal(masterPointi);
1046 
1047  if (!collapsePointToLocation.found(masterPointi))
1048  {
1049  const label otherVertex = e.otherVertex(masterPointi);
1050 
1051  if (!collapsePointToLocation.found(otherVertex))
1052  {
1054  << masterPointi << " on edge " << edgeI << " " << e
1055  << " is not marked for collapse."
1056  << abort(FatalError);
1057  }
1058  else
1059  {
1060  masterPointi = otherVertex;
1061  masterPointPriority = pointPriority[masterPointi];
1062  index = globalPoints.toGlobal(masterPointi);
1063  }
1064  }
1065 
1066  const point& collapsePoint = collapsePointToLocation[masterPointi];
1067 
1068  const pointEdgeCollapse pec
1069  (
1070  collapsePoint,
1071  index,
1072  masterPointPriority
1073  );
1074 
1075  // Mark as collapsible but with nonsense master so it gets
1076  // overwritten and starts an update wave
1077  allEdgeInfo[edgeI] = pointEdgeCollapse
1078  (
1079  collapsePoint,
1080  labelMax,
1081  labelMin
1082  );
1083 
1084  initPointInfo.append(pec);
1085  initPoints.append(e.start());
1086 
1087  initPointInfo.append(pec);
1088  initPoints.append(e.end());
1089 
1090  nCollapsed++;
1091  }
1092  }
1093 
1094  PointEdgeWave<pointEdgeCollapse> collapsePropagator
1095  (
1096  mesh_,
1097  initPoints,
1098  initPointInfo,
1099  allPointInfo,
1100  allEdgeInfo,
1101  mesh_.globalData().nTotalPoints() // Maximum number of iterations
1102  );
1103 
1104  return nCollapsed;
1105 }
1106 
1107 
1108 void Foam::edgeCollapser::filterFace
1109 (
1110  const Map<DynamicList<label>>& collapseStrings,
1111  const List<pointEdgeCollapse>& allPointInfo,
1112  face& f
1113 ) const
1114 {
1115  label newFp = 0;
1116 
1117  face oldFace = f;
1118 
1119  forAll(f, fp)
1120  {
1121  label pointi = f[fp];
1122 
1123  label collapseIndex = allPointInfo[pointi].collapseIndex();
1124 
1125  // Do we have a local point for this index?
1126  if (collapseStrings.found(collapseIndex))
1127  {
1128  label localPointi = collapseStrings[collapseIndex][0];
1129 
1130  if (findIndex(SubList<label>(f, newFp), localPointi) == -1)
1131  {
1132  f[newFp++] = localPointi;
1133  }
1134  }
1135  else if (collapseIndex == -1)
1136  {
1138  << "Point " << pointi << " was not visited by PointEdgeWave"
1139  << endl;
1140  }
1141  else
1142  {
1143  f[newFp++] = pointi;
1144  }
1145  }
1146 
1147 
1148  // Check for pinched face. Tries to correct
1149  // - consecutive duplicate vertex. Removes duplicate vertex.
1150  // - duplicate vertex with one other vertex in between (spike).
1151  // Both of these should not really occur! and should be checked before
1152  // collapsing edges.
1153 
1154  const label size = newFp;
1155 
1156  newFp = 2;
1157 
1158  for (label fp = 2; fp < size; fp++)
1159  {
1160  label fp1 = fp-1;
1161  label fp2 = fp-2;
1162 
1163  label pointi = f[fp];
1164 
1165  // Search for previous occurrence.
1166  label index = findIndex(SubList<label>(f, fp), pointi);
1167 
1168  if (index == fp1)
1169  {
1171  << "Removing consecutive duplicate vertex in face "
1172  << f << endl;
1173  // Don't store current pointi
1174  }
1175  else if (index == fp2)
1176  {
1178  << "Removing non-consecutive duplicate vertex in face "
1179  << f << endl;
1180  // Don't store current pointi and remove previous
1181  newFp--;
1182  }
1183  else if (index != -1)
1184  {
1186  << "Pinched face " << f << endl;
1187  f[newFp++] = pointi;
1188  }
1189  else
1190  {
1191  f[newFp++] = pointi;
1192  }
1193  }
1194 
1195  f.setSize(newFp);
1196 }
1197 
1198 
1199 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1200 
1202 :
1203  mesh_(mesh),
1204  guardFraction_(0),
1205  maxCollapseFaceToPointSideLengthCoeff_(0),
1206  allowEarlyCollapseToPoint_(false),
1207  allowEarlyCollapseCoeff_(0)
1208 {}
1209 
1210 
1212 (
1213  const polyMesh& mesh,
1214  const dictionary& dict
1215 )
1216 :
1217  mesh_(mesh),
1218  guardFraction_
1219  (
1220  dict.lookupOrDefault<scalar>("guardFraction", 0)
1221  ),
1222  maxCollapseFaceToPointSideLengthCoeff_
1223  (
1224  dict.lookupOrDefault<scalar>("maxCollapseFaceToPointSideLengthCoeff", 0)
1225  ),
1226  allowEarlyCollapseToPoint_
1227  (
1228  dict.lookupOrDefault<Switch>("allowEarlyCollapseToPoint", true)
1229  ),
1230  allowEarlyCollapseCoeff_
1231  (
1232  dict.lookupOrDefault<scalar>("allowEarlyCollapseCoeff", 0)
1233  )
1234 {
1235  if (debug)
1236  {
1237  Info<< "Edge Collapser Settings:" << nl
1238  << " Guard Fraction = " << guardFraction_ << nl
1239  << " Max collapse face to point side length = "
1240  << maxCollapseFaceToPointSideLengthCoeff_ << nl
1241  << " " << (allowEarlyCollapseToPoint_ ? "Allow" : "Do not allow")
1242  << " early collapse to point" << nl
1243  << " Early collapse coeff = " << allowEarlyCollapseCoeff_
1244  << endl;
1245  }
1246 }
1247 
1248 
1249 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1250 
1252 (
1253  const List<pointEdgeCollapse>& allPointInfo,
1254  polyTopoChange& meshMod
1255 ) const
1256 {
1257  const cellList& cells = mesh_.cells();
1258  const labelList& faceOwner = mesh_.faceOwner();
1259  const labelList& faceNeighbour = mesh_.faceNeighbour();
1260  const labelListList& pointFaces = mesh_.pointFaces();
1261 
1262 // // Dump point collapses
1263 // label count = 0;
1264 // forAll(allPointInfo, ptI)
1265 // {
1266 // const pointEdgeCollapse& pec = allPointInfo[ptI];
1267 //
1268 // if (mesh_.points()[ptI] != pec.collapsePoint())
1269 // {
1270 // count++;
1271 // }
1272 // }
1273 //
1274 // OFstream str("collapses_" + name(count) + ".obj");
1275 // // Dump point collapses
1276 // forAll(allPointInfo, ptI)
1277 // {
1278 // const pointEdgeCollapse& pec = allPointInfo[ptI];
1279 //
1280 // if
1281 // (
1282 // mesh_.points()[ptI] != pec.collapsePoint()
1283 // && pec.collapsePoint() != vector(great, great, great)
1284 // )
1285 // {
1286 // meshTools::writeOBJ
1287 // (
1288 // str,
1289 // mesh_.points()[ptI],
1290 // pec.collapsePoint()
1291 // );
1292 // }
1293 // }
1294 
1295 
1296 
1297  bool meshChanged = false;
1298 
1299  PackedBoolList removedPoints(mesh_.nPoints());
1300 
1301  // Create strings of edges.
1302  // Map from collapseIndex(=global master point) to set of points
1303  Map<DynamicList<label>> collapseStrings;
1304  {
1305  // 1. Count elements per collapseIndex
1306  Map<label> nPerIndex(mesh_.nPoints()/10);
1307  forAll(allPointInfo, pointi)
1308  {
1309  label collapseIndex = allPointInfo[pointi].collapseIndex();
1310 
1311  if (collapseIndex != -1 && collapseIndex != -2)
1312  {
1313  Map<label>::iterator fnd = nPerIndex.find(collapseIndex);
1314  if (fnd != nPerIndex.end())
1315  {
1316  fnd()++;
1317  }
1318  else
1319  {
1320  nPerIndex.insert(collapseIndex, 1);
1321  }
1322  }
1323  }
1324 
1325  // 2. Size
1326  collapseStrings.resize(2*nPerIndex.size());
1327  forAllConstIter(Map<label>, nPerIndex, iter)
1328  {
1329  collapseStrings.insert(iter.key(), DynamicList<label>(iter()));
1330  }
1331 
1332  // 3. Fill
1333  forAll(allPointInfo, pointi)
1334  {
1335  const label collapseIndex = allPointInfo[pointi].collapseIndex();
1336 
1337  if (collapseIndex != -1 && collapseIndex != -2)
1338  {
1339  collapseStrings[collapseIndex].append(pointi);
1340  }
1341  }
1342  }
1343 
1344 
1345 
1346 
1347 // OFstream str2("collapseStrings_" + name(count) + ".obj");
1348 // // Dump point collapses
1349 // forAllConstIter(Map<DynamicList<label>>, collapseStrings, iter)
1350 // {
1351 // const label masterPoint = iter.key();
1352 // const DynamicList<label>& edgeCollapses = iter();
1353 //
1354 // forAll(edgeCollapses, eI)
1355 // {
1356 // meshTools::writeOBJ
1357 // (
1358 // str2,
1359 // mesh_.points()[edgeCollapses[eI]],
1360 // mesh_.points()[masterPoint]
1361 // );
1362 // }
1363 // }
1364 
1365 
1366 
1367  // Current faces (is also collapseStatus: f.size() < 3)
1368  faceList newFaces(mesh_.faces());
1369 
1370  // Current cellCollapse status
1371  boolList cellRemoved(mesh_.nCells(), false);
1372 
1373  label nUnvisited = 0;
1374  label nUncollapsed = 0;
1375  label nCollapsed = 0;
1376 
1377  forAll(allPointInfo, pI)
1378  {
1379  const pointEdgeCollapse& pec = allPointInfo[pI];
1380 
1381  if (pec.collapseIndex() == -1)
1382  {
1383  nUnvisited++;
1384  }
1385  else if (pec.collapseIndex() == -2)
1386  {
1387  nUncollapsed++;
1388  }
1389  else
1390  {
1391  nCollapsed++;
1392  }
1393  }
1394 
1395  label nPoints = allPointInfo.size();
1396 
1398  reduce(nUnvisited, sumOp<label>());
1399  reduce(nUncollapsed, sumOp<label>());
1400  reduce(nCollapsed, sumOp<label>());
1401 
1402  Info<< incrIndent;
1403  Info<< indent << "Number of points : " << nPoints << nl
1404  << indent << "Not visited : " << nUnvisited << nl
1405  << indent << "Not collapsed : " << nUncollapsed << nl
1406  << indent << "Collapsed : " << nCollapsed << nl
1407  << endl;
1408  Info<< decrIndent;
1409 
1410  do
1411  {
1412  forAll(newFaces, facei)
1413  {
1414  filterFace(collapseStrings, allPointInfo, newFaces[facei]);
1415  }
1416 
1417  // Check if faces to be collapsed cause cells to become collapsed.
1418  label nCellCollapsed = 0;
1419 
1420  forAll(cells, celli)
1421  {
1422  if (!cellRemoved[celli])
1423  {
1424  const cell& cFaces = cells[celli];
1425 
1426  label nFaces = cFaces.size();
1427 
1428  forAll(cFaces, i)
1429  {
1430  label facei = cFaces[i];
1431 
1432  if (newFaces[facei].size() < 3)
1433  {
1434  --nFaces;
1435 
1436  if (nFaces < 4)
1437  {
1438  Pout<< "Cell:" << celli
1439  << " uses faces:" << cFaces
1440  << " of which too many are marked for removal:"
1441  << endl
1442  << " ";
1443 
1444 
1445  forAll(cFaces, j)
1446  {
1447  if (newFaces[cFaces[j]].size() < 3)
1448  {
1449  Pout<< ' '<< cFaces[j];
1450  }
1451  }
1452  Pout<< endl;
1453 
1454  cellRemoved[celli] = true;
1455 
1456  // Collapse all edges of cell to nothing
1457 // collapseEdges(cellEdges[celli]);
1458 
1459  nCellCollapsed++;
1460 
1461  break;
1462  }
1463  }
1464  }
1465  }
1466  }
1467 
1468  reduce(nCellCollapsed, sumOp<label>());
1469  Info<< indent << "Collapsing " << nCellCollapsed << " cells" << endl;
1470 
1471  if (nCellCollapsed == 0)
1472  {
1473  break;
1474  }
1475 
1476  } while (true);
1477 
1478 
1479  // Keep track of faces that have been done already.
1480  boolList doneFace(mesh_.nFaces(), false);
1481 
1482  {
1483  // Mark points used.
1484  boolList usedPoint(mesh_.nPoints(), false);
1485 
1486  forAll(cellRemoved, celli)
1487  {
1488  if (cellRemoved[celli])
1489  {
1490  meshMod.removeCell(celli, -1);
1491  }
1492  }
1493 
1494  // Remove faces
1495  forAll(newFaces, facei)
1496  {
1497  const face& f = newFaces[facei];
1498 
1499  if (f.size() < 3)
1500  {
1501  meshMod.removeFace(facei, -1);
1502  meshChanged = true;
1503 
1504  // Mark face as been done.
1505  doneFace[facei] = true;
1506  }
1507  else
1508  {
1509  // Kept face. Mark vertices
1510  forAll(f, fp)
1511  {
1512  usedPoint[f[fp]] = true;
1513  }
1514  }
1515  }
1516 
1517  // Remove unused vertices that have not been marked for removal already
1518  forAll(usedPoint, pointi)
1519  {
1520  if (!usedPoint[pointi])
1521  {
1522  removedPoints[pointi] = true;
1523  meshMod.removePoint(pointi, -1);
1524  meshChanged = true;
1525  }
1526  }
1527  }
1528 
1529  // Modify the point location of the remaining points
1530  forAll(allPointInfo, pointi)
1531  {
1532  const label collapseIndex = allPointInfo[pointi].collapseIndex();
1533  const point& collapsePoint = allPointInfo[pointi].collapsePoint();
1534 
1535  if
1536  (
1537  removedPoints[pointi] == false
1538  && collapseIndex != -1
1539  && collapseIndex != -2
1540  )
1541  {
1542  meshMod.modifyPoint
1543  (
1544  pointi,
1545  collapsePoint,
1546  false
1547  );
1548  }
1549  }
1550 
1551 
1552  const polyBoundaryMesh& boundaryMesh = mesh_.boundaryMesh();
1553 
1554  // Renumber faces that use points
1555  forAll(allPointInfo, pointi)
1556  {
1557  if (removedPoints[pointi] == true)
1558  {
1559  const labelList& changedFaces = pointFaces[pointi];
1560 
1561  forAll(changedFaces, changedFacei)
1562  {
1563  label facei = changedFaces[changedFacei];
1564 
1565  if (!doneFace[facei])
1566  {
1567  doneFace[facei] = true;
1568 
1569  // Get current connectivity
1570  label own = faceOwner[facei];
1571  label nei = -1;
1572  label patchID = -1;
1573 
1574  if (mesh_.isInternalFace(facei))
1575  {
1576  nei = faceNeighbour[facei];
1577  }
1578  else
1579  {
1580  patchID = boundaryMesh.whichPatch(facei);
1581  }
1582 
1583  meshMod.modifyFace
1584  (
1585  newFaces[facei], // face
1586  facei, // facei to change
1587  own, // owner
1588  nei, // neighbour
1589  false, // flipFaceFlux
1590  patchID // patch
1591  );
1592 
1593  meshChanged = true;
1594  }
1595  }
1596  }
1597  }
1598 
1599  return meshChanged;
1600 }
1601 
1602 
1604 (
1605  const globalIndex& globalPoints,
1606  const labelList& pointPriority,
1607  const Map<point>& collapsePointToLocation,
1609  List<pointEdgeCollapse>& allPointInfo,
1610  const bool allowCellCollapse
1611 ) const
1612 {
1613  // Make sure we don't collapse cells
1614  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1615  const faceList faces = mesh_.faces();
1616  const edgeList& edges = mesh_.edges();
1617  const labelListList& faceEdges = mesh_.faceEdges();
1618  const labelListList& pointEdges = mesh_.pointEdges();
1619  const cellList& cells = mesh_.cells();
1620 
1621  labelHashSet uniqueCollapses;
1622  labelHashSet duplicateCollapses;
1623 
1624  while (true)
1625  {
1626  label nUncollapsed = 0;
1627 
1629  (
1630  mesh_,
1631  collapseEdge,
1633  0
1634  );
1635 
1636  // Create consistent set of collapses
1637  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1638  // Note: requires collapseEdge to be synchronised.
1639  syncCollapse
1640  (
1641  globalPoints,
1642  pointPriority,
1643  collapseEdge,
1644  collapsePointToLocation,
1645  allPointInfo
1646  );
1647 
1648  // Get collapsed faces
1649 
1650  PackedBoolList isCollapsedFace(mesh_.nFaces());
1651  PackedBoolList markedPoints(mesh_.nPoints());
1652 
1653  forAll(faces, facei)
1654  {
1655  const face& f = faces[facei];
1656 
1657  isCollapsedFace[facei] = isFaceCollapsed(f, allPointInfo);
1658 
1659  if (isCollapsedFace[facei] < 1)
1660  {
1661  determineDuplicatePointsOnFace
1662  (
1663  f,
1664  markedPoints,
1665  uniqueCollapses,
1666  duplicateCollapses,
1667  allPointInfo
1668  );
1669  }
1670  }
1671 
1672  // Synchronise the marked points
1674  (
1675  mesh_,
1676  markedPoints,
1678  0
1679  );
1680 
1681  // Mark all edges attached to the point for collapse
1682  forAll(markedPoints, pointi)
1683  {
1684  if (markedPoints[pointi])
1685  {
1686  const label index = allPointInfo[pointi].collapseIndex();
1687 
1688  const labelList& ptEdges = pointEdges[pointi];
1689 
1690  forAll(ptEdges, ptEdgeI)
1691  {
1692  const label edgeI = ptEdges[ptEdgeI];
1693  const label nbrPointi = edges[edgeI].otherVertex(pointi);
1694  const label nbrIndex
1695  = allPointInfo[nbrPointi].collapseIndex();
1696 
1697  if (collapseEdge[edgeI] && nbrIndex == index)
1698  {
1699  collapseEdge[edgeI] = false;
1700  nUncollapsed++;
1701  }
1702  }
1703  }
1704  }
1705 
1706  if (!allowCellCollapse)
1707  {
1708  // Check collapsed cells
1709  forAll(cells, celli)
1710  {
1711  const cell& cFaces = cells[celli];
1712 
1713  label nFaces = cFaces.size();
1714 
1715  forAll(cFaces, fI)
1716  {
1717  label facei = cFaces[fI];
1718 
1719  if (isCollapsedFace[facei])
1720  {
1721  nFaces--;
1722  }
1723  }
1724 
1725  if (nFaces < 4)
1726  {
1727  forAll(cFaces, fI)
1728  {
1729  label facei = cFaces[fI];
1730 
1731  const labelList& fEdges = faceEdges[facei];
1732 
1733  // Unmark this face for collapse
1734  forAll(fEdges, fEdgeI)
1735  {
1736  label edgeI = fEdges[fEdgeI];
1737 
1738  if (collapseEdge[edgeI])
1739  {
1740  collapseEdge[edgeI] = false;
1741  nUncollapsed++;
1742  }
1743  }
1744 
1745  nFaces += isCollapsedFace[facei] ? 1 : 0;
1746 
1747  // Uncollapsed this face.
1748  isCollapsedFace[facei] = false;
1749  }
1750  }
1751 
1752  if (nFaces < 4)
1753  {
1755  << "Cell " << celli << " " << cFaces << nl
1756  << "is " << nFaces << ", "
1757  << "but cell collapse has been disabled."
1758  << abort(FatalError);
1759  }
1760  }
1761  }
1762 
1763  nUncollapsed += breakStringsAtEdges(collapseEdge, allPointInfo);
1764 
1765  reduce(nUncollapsed, sumOp<label>());
1766 
1767  Info<< " Uncollapsed edges = " << nUncollapsed << " / "
1768  << returnReduce(mesh_.nEdges(), sumOp<label>()) << endl;
1769 
1770  if (nUncollapsed == 0)
1771  {
1772  break;
1773  }
1774  }
1775 }
1776 
1777 
1779 (
1780  const scalarField& minEdgeLen,
1781  const labelList& pointPriority,
1783  Map<point>& collapsePointToLocation
1784 ) const
1785 {
1786  // Work out which edges to collapse
1787  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1788 
1789  const pointField& points = mesh_.points();
1790  const edgeList& edges = mesh_.edges();
1791 
1792  label nCollapsed = 0;
1793 
1794  forAll(edges, edgeI)
1795  {
1796  const edge& e = edges[edgeI];
1797 
1798  if (!collapseEdge[edgeI])
1799  {
1800  if (e.mag(points) < minEdgeLen[edgeI])
1801  {
1802  collapseEdge[edgeI] = true;
1803 
1804  label masterPointi = edgeMaster(pointPriority, e);
1805 
1806  if (masterPointi == -1)
1807  {
1808  const point average
1809  = 0.5*(points[e.start()] + points[e.end()]);
1810 
1811  collapsePointToLocation.set(e.start(), average);
1812  }
1813  else
1814  {
1815  const point& collapsePt = points[masterPointi];
1816 
1817  collapsePointToLocation.set(masterPointi, collapsePt);
1818  }
1819 
1820 
1821  nCollapsed++;
1822  }
1823  }
1824  }
1825 
1826  return nCollapsed;
1827 }
1828 
1829 
1831 (
1832  const scalar maxCos,
1833  const labelList& pointPriority,
1835  Map<point>& collapsePointToLocation
1836 ) const
1837 {
1838  const edgeList& edges = mesh_.edges();
1839  const pointField& points = mesh_.points();
1840  const labelListList& pointEdges = mesh_.pointEdges();
1841 
1842  // Point removal engine
1843  removePoints pointRemover(mesh_, false);
1844 
1845  // Find out points that can be deleted
1846  boolList pointCanBeDeleted;
1847  label nTotRemove = pointRemover.countPointUsage(maxCos, pointCanBeDeleted);
1848 
1849  // Rework point-to-remove into edge-to-collapse.
1850 
1851  label nCollapsed = 0;
1852 
1853  if (nTotRemove > 0)
1854  {
1855  forAll(pointEdges, pointi)
1856  {
1857  if (pointCanBeDeleted[pointi])
1858  {
1859  const labelList& pEdges = pointEdges[pointi];
1860 
1861  if (pEdges.size() == 2)
1862  {
1863  // Always the case?
1864 
1865  label e0 = pEdges[0];
1866  label e1 = pEdges[1];
1867 
1868  if (!collapseEdge[e0] && !collapseEdge[e1])
1869  {
1870  // Get lengths of both edges and choose the smallest
1871  scalar e0length = mag
1872  (
1873  points[edges[e0][0]] - points[edges[e0][1]]
1874  );
1875 
1876  scalar e1length = mag
1877  (
1878  points[edges[e1][0]] - points[edges[e1][1]]
1879  );
1880 
1881  if (e0length <= e1length)
1882  {
1883  collapseEdge[e0] = true;
1884 
1885  checkBoundaryPointMergeEdges
1886  (
1887  pointi,
1888  edges[e0].otherVertex(pointi),
1889  pointPriority,
1890  collapsePointToLocation
1891  );
1892  }
1893  else
1894  {
1895  collapseEdge[e1] = true;
1896 
1897  checkBoundaryPointMergeEdges
1898  (
1899  pointi,
1900  edges[e1].otherVertex(pointi),
1901  pointPriority,
1902  collapsePointToLocation
1903  );
1904  }
1905 
1906  nCollapsed++;
1907  }
1908  }
1909  }
1910  }
1911  }
1912 
1913  return nCollapsed;
1914 }
1915 
1916 
1918 (
1919  const scalarField& faceFilterFactor,
1920  const labelList& pointPriority,
1922  Map<point>& collapsePointToLocation
1923 ) const
1924 {
1925  const faceList& faces = mesh_.faces();
1926 
1927  const scalarField targetFaceSizes = calcTargetFaceSizes();
1928 
1929  // Calculate number of faces that will be collapsed to a point or an edge
1930  label nCollapseToPoint = 0;
1931  label nCollapseToEdge = 0;
1932 
1933  forAll(faces, fI)
1934  {
1935  const face& f = faces[fI];
1936 
1937  if (faceFilterFactor[fI] <= 0)
1938  {
1939  continue;
1940  }
1941 
1942  collapseType flagCollapseFace = collapseFace
1943  (
1944  pointPriority,
1945  f,
1946  fI,
1947  targetFaceSizes[fI],
1948  collapseEdge,
1949  collapsePointToLocation,
1950  faceFilterFactor
1951  );
1952 
1953  if (flagCollapseFace == noCollapse)
1954  {
1955  continue;
1956  }
1957  else if (flagCollapseFace == toPoint)
1958  {
1959  nCollapseToPoint++;
1960  }
1961  else if (flagCollapseFace == toEdge)
1962  {
1963  nCollapseToEdge++;
1964  }
1965  else
1966  {
1968  << "Face is marked to be collapsed to " << flagCollapseFace
1969  << ". Currently can only collapse to point/edge."
1970  << abort(FatalError);
1971  }
1972  }
1973 
1974  return labelPair(nCollapseToPoint, nCollapseToEdge);
1975 }
1976 
1977 
1979 (
1980  const faceZone& fZone,
1981  const scalarField& faceFilterFactor,
1982  const labelList& pointPriority,
1984  Map<point>& collapsePointToLocation
1985 ) const
1986 {
1987  const faceList& faces = mesh_.faces();
1988 
1989  const scalarField targetFaceSizes = calcTargetFaceSizes();
1990 
1991  // Calculate number of faces that will be collapsed to a point or an edge
1992  label nCollapseToPoint = 0;
1993  label nCollapseToEdge = 0;
1994 
1995  forAll(faces, fI)
1996  {
1997  if (fZone.localIndex(fI) == -1)
1998  {
1999  continue;
2000  }
2001 
2002  const face& f = faces[fI];
2003 
2004  if (faceFilterFactor[fI] <= 0)
2005  {
2006  continue;
2007  }
2008 
2009  collapseType flagCollapseFace = collapseFace
2010  (
2011  pointPriority,
2012  f,
2013  fI,
2014  targetFaceSizes[fI],
2015  collapseEdge,
2016  collapsePointToLocation,
2017  faceFilterFactor
2018  );
2019 
2020  if (flagCollapseFace == noCollapse)
2021  {
2022  continue;
2023  }
2024  else if (flagCollapseFace == toPoint)
2025  {
2026  nCollapseToPoint++;
2027  }
2028  else if (flagCollapseFace == toEdge)
2029  {
2030  nCollapseToEdge++;
2031  }
2032  else
2033  {
2035  << "Face is marked to be collapsed to " << flagCollapseFace
2036  << ". Currently can only collapse to point/edge."
2037  << abort(FatalError);
2038  }
2039  }
2040 
2041  return labelPair(nCollapseToPoint, nCollapseToEdge);
2042 
2043 // const edgeList& edges = mesh_.edges();
2044 // const pointField& points = mesh_.points();
2045 // const labelListList& edgeFaces = mesh_.edgeFaces();
2046 // const polyBoundaryMesh& bMesh = mesh_.boundaryMesh();
2047 //
2048 // forAll(edges, eI)
2049 // {
2050 // const edge& e = edges[eI];
2051 //
2052 // const labelList& eFaces = edgeFaces[eI];
2053 //
2054 // bool keepEdge = false;
2055 //
2056 // label nInternalFaces = 0;
2057 // label nPatchFaces = 0;
2058 // label nIndirectFaces = 0;
2059 //
2060 // bool coupled = false;
2061 //
2062 // forAll(eFaces, eFacei)
2063 // {
2064 // const label eFaceIndex = eFaces[eFacei];
2065 //
2066 // if (mesh_.isInternalFace(eFaceIndex))
2067 // {
2068 // nInternalFaces++;
2069 // }
2070 // else
2071 // {
2072 // const label patchIndex = bMesh.whichPatch(eFaceIndex);
2073 // const polyPatch& pPatch = bMesh[patchIndex];
2074 //
2075 // if (pPatch.coupled())
2076 // {
2077 // coupled = true;
2078 // nInternalFaces++;
2079 // }
2080 // else
2081 // {
2082 // // Keep the edge if an attached face is not in the zone
2083 // if (fZone.localIndex(eFaceIndex) == -1)
2084 // {
2085 // nPatchFaces++;
2086 // }
2087 // else
2088 // {
2089 // nIndirectFaces++;
2090 // }
2091 // }
2092 // }
2093 // }
2094 //
2095 // if (eFaces.size() != nInternalFaces + nPatchFaces + nIndirectFaces)
2096 // {
2097 // Pout<< eFaces.size() << " ("
2098 // << nInternalFaces << "/" << nPatchFaces << "/"
2099 // << nIndirectFaces << ")" << endl;
2100 // }
2101 //
2102 // if
2103 // (
2104 // eFaces.size() == nInternalFaces
2105 // || nIndirectFaces < (coupled ? 1 : 2)
2106 // )
2107 // {
2108 // keepEdge = true;
2109 // }
2110 //
2111 // if (!keepEdge)
2112 // {
2113 // collapseEdge[eI] = true;
2114 //
2115 // const Foam::point collapsePoint =
2116 // 0.5*(points[e.end()] + points[e.start()]);
2117 //
2118 // collapsePointToLocation.insert(e.start(), collapsePoint);
2119 // collapsePointToLocation.insert(e.end(), collapsePoint);
2120 // }
2121 // }
2122 
2123 // OFstream str
2124 // (
2125 // mesh_.time().path()
2126 // /"markedEdges_" + name(collapseEdge.count()) + ".obj"
2127 // );
2128 // label count = 0;
2129 //
2130 // forAll(collapseEdge, eI)
2131 // {
2132 // if (collapseEdge[eI])
2133 // {
2134 // const edge& e = edges[eI];
2135 //
2136 // meshTools::writeOBJ
2137 // (
2138 // str,
2139 // points[e.start()],
2140 // points[e.end()],
2141 // count
2142 // );
2143 // }
2144 // }
2145 }
2146 
2147 
2148 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
bool set(const Key &, const T &newElmt)
Set a new hashedEntry, overwriting existing entries.
Definition: HashTableI.H:91
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
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:167
void resize(const label newSize)
Resize the hash table for efficiency.
Definition: HashTable.C:506
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
A HashTable to objects of type <T> with a label key.
Definition: Map.H:52
A bit-packed bool list.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
static direction size()
Return the number of elements in the VectorSpace = Ncmpts.
Definition: VectorSpaceI.H:83
label localIndex(const label globalIndex) const
Map storing the local index for every global index. Used to find.
Definition: Zone.C:240
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Does polyTopoChanges to remove edges. Can remove faces due to edge collapse but can not remove cells ...
Definition: edgeCollapser.H:67
static label checkMeshQuality(const polyMesh &mesh, const dictionary &meshQualityDict, PackedBoolList &isErrorPoint)
Check mesh and mark points on faces in error.
Definition: edgeCollapser.C:82
labelPair markSmallSliverFaces(const scalarField &faceFilterFactor, const labelList &pointPriority, PackedBoolList &collapseEdge, Map< point > &collapsePointToLocation) const
Find small faces and sliver faces in the mesh and mark the.
edgeCollapser(const polyMesh &mesh)
Construct from mesh.
label markSmallEdges(const scalarField &minEdgeLen, const labelList &pointPriority, PackedBoolList &collapseEdge, Map< point > &collapsePointToLocation) const
Mark (in collapseEdge) any edges to collapse.
labelPair markFaceZoneEdges(const faceZone &fZone, const scalarField &faceFilterFactor, const labelList &pointPriority, PackedBoolList &collapseEdge, Map< point > &collapsePointToLocation) const
Marks edges in the faceZone indirectPatchFaces for collapse.
bool setRefinement(const List< pointEdgeCollapse > &allPointInfo, polyTopoChange &meshMod) const
Play commands into polyTopoChange to create mesh.
void consistentCollapse(const globalIndex &globalPoints, const labelList &pointPriority, const Map< point > &collapsePointToLocation, PackedBoolList &collapseEdge, List< pointEdgeCollapse > &allPointInfo, const bool allowCellCollapse=false) const
Ensure that the collapse is parallel consistent and update.
label markMergeEdges(const scalar maxCos, const labelList &pointPriority, PackedBoolList &collapseEdge, Map< point > &collapsePointToLocation) const
Mark (in collapseEdge) any edges to merge.
static labelHashSet checkBadFaces(const polyMesh &mesh, const dictionary &meshQualityDict)
Calls meshCheck::checkMesh and returns a set of bad faces.
Definition: edgeCollapser.C:46
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:101
Determines length of string of edges walked to point.
Foam::polyBoundaryMesh.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
Direct mesh changes based on v1.3 polyTopoChange syntax.
void modifyPoint(const label, const point &, const bool inCell)
Modify coordinate.
void removePoint(const label, const label)
Remove point / merge points.
void removeFace(const label, const label)
Remove face / merge faces.
void removeCell(const label, const label)
Remove cell / merge cells.
void modifyFace(const face &f, const label facei, const label own, const label nei, const bool flipFaceFlux, const label patchID)
Modify vertices or cell of face.
label nFaces() const
const vectorField & faceAreas() const
Removes selected points from mesh and updates faces using these points.
Definition: removePoints.H:59
label countPointUsage(const scalar minCos, boolList &pointCanBeDeleted) const
Mark in pointCanBeDeleted the points that can be deleted.
Definition: removePoints.C:125
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:436
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
label collapseEdge(triSurface &surf, const scalar minLen)
Keep collapsing all edges < minLen.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const pointField & points
label nPoints
const cellShapeList & cells
const fvPatchList & patches
Functions for checking mesh topology and geometry.
#define WarningInFunction
Report a warning using Foam::Warning.
bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet &wrongFaces)
Check (subset of mesh including baffles) with mesh settings in dict.
Definition: checkMesh.C:33
Namespace for OpenFOAM.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
static const zero Zero
Definition: zero.H:97
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
const doubleScalar e
Definition: doubleScalar.H:106
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensionedTensor eigenVectors(const dimensionedTensor &dt)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
vector point
Point is a vector.
Definition: point.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedVector eigenValues(const dimensionedTensor &dt)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
label longestEdge(const face &f, const pointField &pts)
Find the longest edge on a face. Face point labels index into pts.
Definition: face.C:467
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
List< edge > edgeList
Definition: edgeList.H:38
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
static const label labelMax
Definition: label.H:62
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
UList< label > labelUList
Definition: UList.H:65
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const label labelMin
Definition: label.H:61
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
label newPointi
Definition: readKivaGrid.H:495
labelList f(nPoints)
labelList pointLabels(nPoints, -1)
dictionary dict
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112