edgeCollapser.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-2014 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 "motionSmoother.H"
35 #include "OFstream.H"
36 
37 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 defineTypeNameAndDebug(edgeCollapser, 0);
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(vector::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(vector::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(vector::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 = eigenVector(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
501  forAll(patches, patchI)
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 
534  forAll(patches, patchI)
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 = vector::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 // SeriousErrorIn("collapseFace")
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  {
665  WarningIn
666  (
667  "Foam::conformalVoronoiMesh::collapseFace"
668  )
669  << "All points on one side of face centre, not collapsing."
670  << endl;
671  }
672 
673 // Info<< "Face : " << f << nl
674 // << " Collapse Axis: " << collapseAxis << nl
675 // << " Aspect Ratio : " << aspectRatio << endl;
676 
677  collapseType typeOfCollapse = noCollapse;
678 
679  if (magSqr(collapseAxis) < VSMALL)
680  {
681  typeOfCollapse = toPoint;
682  }
683  else if (fA < aspectRatio*sqr(targetFaceSize*collapseSizeLimitCoeff))
684  {
685  if
686  (
687  allowEarlyCollapseToPoint_
688  && (d.last() - d.first())
689  < targetFaceSize
690  *allowEarlyCollapseCoeff_*maxCollapseFaceToPointSideLengthCoeff_
691  )
692  {
693  typeOfCollapse = toPoint;
694  }
695  else if
696  (
697  (dNeg.last() < guardFraction_*dNeg.first())
698  && (dPos.first() > guardFraction_*dPos.last())
699  )
700  {
701  typeOfCollapse = toEdge;
702  }
703  else if
704  (
705  (d.last() - d.first())
706  < targetFaceSize
707  *maxCollapseFaceToPointSideLengthCoeff_
708  )
709  {
710  // If the face can't be collapsed to an edge, and it has a
711  // small enough span, collapse it to a point.
712  typeOfCollapse = toPoint;
713  }
714  }
715 
716  if (typeOfCollapse == toPoint)
717  {
718  collapseToPoint
719  (
720  faceI,
721  pts,
722  pointPriority,
723  fC,
724  facePts,
725  collapseEdge,
726  collapsePointToLocation
727  );
728  }
729  else if (typeOfCollapse == toEdge)
730  {
731  collapseToEdge
732  (
733  faceI,
734  pts,
735  pointPriority,
736  collapseAxis,
737  fC,
738  facePtsNeg,
739  facePtsPos,
740  dNeg,
741  dPos,
742  dShift,
743  collapseEdge,
744  collapsePointToLocation
745  );
746  }
747 
748  return typeOfCollapse;
749 }
750 
751 
752 Foam::label Foam::edgeCollapser::edgeMaster
753 (
754  const labelList& pointPriority,
755  const edge& e
756 ) const
757 {
758  label masterPoint = -1;
759 
760  const label e0 = e.start();
761  const label e1 = e.end();
762 
763  const label e0Priority = pointPriority[e0];
764  const label e1Priority = pointPriority[e1];
765 
766  if (e0Priority > e1Priority)
767  {
768  masterPoint = e0;
769  }
770  else if (e0Priority < e1Priority)
771  {
772  masterPoint = e1;
773  }
774  else if (e0Priority == e1Priority)
775  {
776  masterPoint = e0;
777  }
778 
779 // // Collapse edge to point with higher priority.
780 // if (pointPriority[e0] >= 0)
781 // {
782 // if (pointPriority[e1] >= 0)
783 // {
784 // // Both points have high priority. Choose one to collapse to.
785 // // Note: should look at feature edges/points!
786 // masterPoint = e0;
787 // }
788 // else
789 // {
790 // masterPoint = e0;
791 // }
792 // }
793 // else
794 // {
795 // if (pointPriority[e1] >= 0)
796 // {
797 // masterPoint = e1;
798 // }
799 // else
800 // {
801 // // None on boundary. Neither is a master.
802 // return -1;
803 // }
804 // }
805 
806  return masterPoint;
807 }
808 
809 
810 void Foam::edgeCollapser::checkBoundaryPointMergeEdges
811 (
812  const label pointI,
813  const label otherPointI,
814  const labelList& pointPriority,
815  Map<point>& collapsePointToLocation
816 ) const
817 {
818  const pointField& points = mesh_.points();
819 
820  const label e0Priority = pointPriority[pointI];
821  const label e1Priority = pointPriority[otherPointI];
822 
823  if (e0Priority > e1Priority)
824  {
825  collapsePointToLocation.set
826  (
827  otherPointI,
828  points[pointI]
829  );
830  }
831  else if (e0Priority < e1Priority)
832  {
833  collapsePointToLocation.set
834  (
835  pointI,
836  points[otherPointI]
837  );
838  }
839  else // e0Priority == e1Priority
840  {
841  collapsePointToLocation.set
842  (
843  pointI,
844  points[otherPointI]
845  );
846 
847 // Foam::point averagePt
848 // (
849 // 0.5*(points[otherPointI] + points[pointI])
850 // );
851 //
852 // collapsePointToLocation.set(pointI, averagePt);
853 // collapsePointToLocation.set(otherPointI, averagePt);
854  }
855 }
856 
857 
858 Foam::label Foam::edgeCollapser::breakStringsAtEdges
859 (
860  const PackedBoolList& markedEdges,
861  PackedBoolList& collapseEdge,
862  List<pointEdgeCollapse>& allPointInfo
863 ) const
864 {
865  const edgeList& edges = mesh_.edges();
866  const labelListList& pointEdges = mesh_.pointEdges();
867 
868  label nUncollapsed = 0;
869 
870  forAll(edges, eI)
871  {
872  if (markedEdges[eI])
873  {
874  const edge& e = edges[eI];
875 
876  const label startCollapseIndex
877  = allPointInfo[e.start()].collapseIndex();
878 
879  if (startCollapseIndex != -1 && startCollapseIndex != -2)
880  {
881  const label endCollapseIndex
882  = allPointInfo[e.end()].collapseIndex();
883 
884  if
885  (
886  !collapseEdge[eI]
887  && startCollapseIndex == endCollapseIndex
888  )
889  {
890  const labelList& ptEdgesStart = pointEdges[e.start()];
891 
892  forAll(ptEdgesStart, ptEdgeI)
893  {
894  const label edgeI = ptEdgesStart[ptEdgeI];
895 
896  const label nbrPointI
897  = edges[edgeI].otherVertex(e.start());
898  const label nbrIndex
899  = allPointInfo[nbrPointI].collapseIndex();
900 
901  if
902  (
903  collapseEdge[edgeI]
904  && nbrIndex == startCollapseIndex
905  )
906  {
907  collapseEdge[edgeI] = false;
908  nUncollapsed++;
909  }
910  }
911  }
912  }
913  }
914  }
915 
916  return nUncollapsed;
917 }
918 
919 
920 void Foam::edgeCollapser::determineDuplicatePointsOnFace
921 (
922  const face& f,
923  PackedBoolList& markedPoints,
924  labelHashSet& uniqueCollapses,
925  labelHashSet& duplicateCollapses,
926  List<pointEdgeCollapse>& allPointInfo
927 ) const
928 {
929  uniqueCollapses.clear();
930  duplicateCollapses.clear();
931 
932  forAll(f, fpI)
933  {
934  label index = allPointInfo[f[fpI]].collapseIndex();
935 
936  // Check for consecutive duplicate
937  if (index != allPointInfo[f.prevLabel(fpI)].collapseIndex())
938  {
939  if (!uniqueCollapses.insert(index))
940  {
941  // Failed inserting so duplicate
942  duplicateCollapses.insert(index);
943  }
944  }
945  }
946 
947  // Now duplicateCollapses contains duplicate collapse indices.
948  // Convert to points.
949  forAll(f, fpI)
950  {
951  label index = allPointInfo[f[fpI]].collapseIndex();
952  if (duplicateCollapses.found(index))
953  {
954  markedPoints[f[fpI]] = true;
955  }
956  }
957 }
958 
959 
960 Foam::label Foam::edgeCollapser::countEdgesOnFace
961 (
962  const face& f,
963  List<pointEdgeCollapse>& allPointInfo
964 ) const
965 {
966  label nEdges = 0;
967 
968  forAll(f, fpI)
969  {
970  const label pointI = f[fpI];
971  const label newPointI = allPointInfo[pointI].collapseIndex();
972 
973  if (newPointI == -2)
974  {
975  nEdges++;
976  }
977  else
978  {
979  const label prevPointI = f[f.fcIndex(fpI)];
980  const label prevNewPointI
981  = allPointInfo[prevPointI].collapseIndex();
982 
983  if (newPointI != prevNewPointI)
984  {
985  nEdges++;
986  }
987  }
988  }
989 
990  return nEdges;
991 }
992 
993 
994 bool Foam::edgeCollapser::isFaceCollapsed
995 (
996  const face& f,
997  List<pointEdgeCollapse>& allPointInfo
998 ) const
999 {
1000  label nEdges = countEdgesOnFace(f, allPointInfo);
1001 
1002  // Polygons must have 3 or more edges to be valid
1003  if (nEdges < 3)
1004  {
1005  return true;
1006  }
1007 
1008  return false;
1009 }
1010 
1011 
1012 // Create consistent set of collapses.
1013 // collapseEdge : per edge:
1014 // -1 : do not collapse
1015 // 0 : collapse to start
1016 // 1 : collapse to end
1017 // Note: collapseEdge has to be parallel consistent (in orientation)
1018 Foam::label Foam::edgeCollapser::syncCollapse
1019 (
1020  const globalIndex& globalPoints,
1021  const labelList& pointPriority,
1022  const PackedBoolList& collapseEdge,
1023  const Map<point>& collapsePointToLocation,
1024  List<pointEdgeCollapse>& allPointInfo
1025 ) const
1026 {
1027  const edgeList& edges = mesh_.edges();
1028 
1029  label nCollapsed = 0;
1030 
1031  DynamicList<label> initPoints(mesh_.nPoints());
1032  DynamicList<pointEdgeCollapse> initPointInfo(mesh_.nPoints());
1033 
1034  allPointInfo.clear();
1035  allPointInfo.setSize(mesh_.nPoints());
1036 
1037  // Initialise edges to no collapse
1038  List<pointEdgeCollapse> allEdgeInfo
1039  (
1040  mesh_.nEdges(),
1042  );
1043 
1044  // Mark selected edges for collapse
1045  forAll(edges, edgeI)
1046  {
1047  if (collapseEdge[edgeI])
1048  {
1049  const edge& e = edges[edgeI];
1050 
1051  label masterPointI = e.start();
1052 
1053  // Choose the point on the edge with the highest priority.
1054  if (pointPriority[e.end()] > pointPriority[e.start()])
1055  {
1056  masterPointI = e.end();
1057  }
1058 
1059  label masterPointPriority = pointPriority[masterPointI];
1060 
1061  label index = globalPoints.toGlobal(masterPointI);
1062 
1063  if (!collapsePointToLocation.found(masterPointI))
1064  {
1065  const label otherVertex = e.otherVertex(masterPointI);
1066 
1067  if (!collapsePointToLocation.found(otherVertex))
1068  {
1069  FatalErrorIn
1070  (
1071  "syncCollapse\n"
1072  "(\n"
1073  " const polyMesh&,\n"
1074  " const globalIndex&,\n"
1075  " const labelList&,\n"
1076  " const PackedBoolList&,\n"
1077  " Map<point>&,\n"
1078  " List<pointEdgeCollapse>&\n"
1079  ")\n"
1080  ) << masterPointI << " on edge " << edgeI << " " << e
1081  << " is not marked for collapse."
1082  << abort(FatalError);
1083  }
1084  else
1085  {
1086  masterPointI = otherVertex;
1087  masterPointPriority = pointPriority[masterPointI];
1088  index = globalPoints.toGlobal(masterPointI);
1089  }
1090  }
1091 
1092  const point& collapsePoint = collapsePointToLocation[masterPointI];
1093 
1094  const pointEdgeCollapse pec
1095  (
1096  collapsePoint,
1097  index,
1098  masterPointPriority
1099  );
1100 
1101  // Mark as collapsable but with nonsense master so it gets
1102  // overwritten and starts an update wave
1103  allEdgeInfo[edgeI] = pointEdgeCollapse
1104  (
1105  collapsePoint,
1106  labelMax,
1107  labelMin
1108  );
1109 
1110  initPointInfo.append(pec);
1111  initPoints.append(e.start());
1112 
1113  initPointInfo.append(pec);
1114  initPoints.append(e.end());
1115 
1116  nCollapsed++;
1117  }
1118  }
1119 
1120  PointEdgeWave<pointEdgeCollapse> collapsePropagator
1121  (
1122  mesh_,
1123  initPoints,
1124  initPointInfo,
1125  allPointInfo,
1126  allEdgeInfo,
1127  mesh_.globalData().nTotalPoints() // Maximum number of iterations
1128  );
1129 
1130  return nCollapsed;
1131 }
1132 
1133 
1134 void Foam::edgeCollapser::filterFace
1135 (
1136  const Map<DynamicList<label> >& collapseStrings,
1137  const List<pointEdgeCollapse>& allPointInfo,
1138  face& f
1139 ) const
1140 {
1141  label newFp = 0;
1142 
1143  face oldFace = f;
1144 
1145  forAll(f, fp)
1146  {
1147  label pointI = f[fp];
1148 
1149  label collapseIndex = allPointInfo[pointI].collapseIndex();
1150 
1151  // Do we have a local point for this index?
1152  if (collapseStrings.found(collapseIndex))
1153  {
1154  label localPointI = collapseStrings[collapseIndex][0];
1155 
1156  if (findIndex(SubList<label>(f, newFp), localPointI) == -1)
1157  {
1158  f[newFp++] = localPointI;
1159  }
1160  }
1161  else if (collapseIndex == -1)
1162  {
1163  WarningIn
1164  (
1165  "filterFace"
1166  "(const label, const Map<DynamicList<label> >&, face&)"
1167  ) << "Point " << pointI << " was not visited by PointEdgeWave"
1168  << endl;
1169  }
1170  else
1171  {
1172  f[newFp++] = pointI;
1173  }
1174  }
1175 
1176 
1177  // Check for pinched face. Tries to correct
1178  // - consecutive duplicate vertex. Removes duplicate vertex.
1179  // - duplicate vertex with one other vertex in between (spike).
1180  // Both of these should not really occur! and should be checked before
1181  // collapsing edges.
1182 
1183  const label size = newFp;
1184 
1185  newFp = 2;
1186 
1187  for (label fp = 2; fp < size; fp++)
1188  {
1189  label fp1 = fp-1;
1190  label fp2 = fp-2;
1191 
1192  label pointI = f[fp];
1193 
1194  // Search for previous occurrence.
1195  label index = findIndex(SubList<label>(f, fp), pointI);
1196 
1197  if (index == fp1)
1198  {
1199  WarningIn
1200  (
1201  "Foam::edgeCollapser::filterFace(const label faceI, "
1202  "face& f) const"
1203  ) << "Removing consecutive duplicate vertex in face "
1204  << f << endl;
1205  // Don't store current pointI
1206  }
1207  else if (index == fp2)
1208  {
1209  WarningIn
1210  (
1211  "Foam::edgeCollapser::filterFace(const label faceI, "
1212  "face& f) const"
1213  ) << "Removing non-consecutive duplicate vertex in face "
1214  << f << endl;
1215  // Don't store current pointI and remove previous
1216  newFp--;
1217  }
1218  else if (index != -1)
1219  {
1220  WarningIn
1221  (
1222  "Foam::edgeCollapser::filterFace(const label faceI, "
1223  "face& f) const"
1224  ) << "Pinched face " << f << endl;
1225  f[newFp++] = pointI;
1226  }
1227  else
1228  {
1229  f[newFp++] = pointI;
1230  }
1231  }
1232 
1233  f.setSize(newFp);
1234 }
1235 
1236 
1237 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1238 
1239 Foam::edgeCollapser::edgeCollapser(const polyMesh& mesh)
1240 :
1241  mesh_(mesh),
1242  guardFraction_(0),
1243  maxCollapseFaceToPointSideLengthCoeff_(0),
1244  allowEarlyCollapseToPoint_(false),
1245  allowEarlyCollapseCoeff_(0)
1246 {}
1247 
1248 
1249 Foam::edgeCollapser::edgeCollapser
1251  const polyMesh& mesh,
1252  const dictionary& dict
1253 )
1254 :
1255  mesh_(mesh),
1256  guardFraction_
1257  (
1258  dict.lookupOrDefault<scalar>("guardFraction", 0)
1259  ),
1260  maxCollapseFaceToPointSideLengthCoeff_
1261  (
1262  dict.lookupOrDefault<scalar>("maxCollapseFaceToPointSideLengthCoeff", 0)
1263  ),
1264  allowEarlyCollapseToPoint_
1265  (
1266  dict.lookupOrDefault<Switch>("allowEarlyCollapseToPoint", true)
1267  ),
1268  allowEarlyCollapseCoeff_
1269  (
1270  dict.lookupOrDefault<scalar>("allowEarlyCollapseCoeff", 0)
1271  )
1272 {
1273  if (debug)
1274  {
1275  Info<< "Edge Collapser Settings:" << nl
1276  << " Guard Fraction = " << guardFraction_ << nl
1277  << " Max collapse face to point side length = "
1278  << maxCollapseFaceToPointSideLengthCoeff_ << nl
1279  << " " << (allowEarlyCollapseToPoint_ ? "Allow" : "Do not allow")
1280  << " early collapse to point" << nl
1281  << " Early collapse coeff = " << allowEarlyCollapseCoeff_
1282  << endl;
1283  }
1284 }
1285 
1286 
1287 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1288 
1291  const List<pointEdgeCollapse>& allPointInfo,
1292  polyTopoChange& meshMod
1293 ) const
1294 {
1295  const cellList& cells = mesh_.cells();
1296  const labelList& faceOwner = mesh_.faceOwner();
1297  const labelList& faceNeighbour = mesh_.faceNeighbour();
1298  const labelListList& pointFaces = mesh_.pointFaces();
1299  const pointZoneMesh& pointZones = mesh_.pointZones();
1300 
1301 
1302 
1303 
1304 // // Dump point collapses
1305 // label count = 0;
1306 // forAll(allPointInfo, ptI)
1307 // {
1308 // const pointEdgeCollapse& pec = allPointInfo[ptI];
1309 //
1310 // if (mesh_.points()[ptI] != pec.collapsePoint())
1311 // {
1312 // count++;
1313 // }
1314 // }
1315 //
1316 // OFstream str("collapses_" + name(count) + ".obj");
1317 // // Dump point collapses
1318 // forAll(allPointInfo, ptI)
1319 // {
1320 // const pointEdgeCollapse& pec = allPointInfo[ptI];
1321 //
1322 // if
1323 // (
1324 // mesh_.points()[ptI] != pec.collapsePoint()
1325 // && pec.collapsePoint() != vector(GREAT, GREAT, GREAT)
1326 // )
1327 // {
1328 // meshTools::writeOBJ
1329 // (
1330 // str,
1331 // mesh_.points()[ptI],
1332 // pec.collapsePoint()
1333 // );
1334 // }
1335 // }
1336 
1337 
1338 
1339  bool meshChanged = false;
1340 
1341  PackedBoolList removedPoints(mesh_.nPoints());
1342 
1343  // Create strings of edges.
1344  // Map from collapseIndex(=global master point) to set of points
1345  Map<DynamicList<label> > collapseStrings;
1346  {
1347  // 1. Count elements per collapseIndex
1348  Map<label> nPerIndex(mesh_.nPoints()/10);
1349  forAll(allPointInfo, pointI)
1350  {
1351  label collapseIndex = allPointInfo[pointI].collapseIndex();
1352 
1353  if (collapseIndex != -1 && collapseIndex != -2)
1354  {
1355  Map<label>::iterator fnd = nPerIndex.find(collapseIndex);
1356  if (fnd != nPerIndex.end())
1357  {
1358  fnd()++;
1359  }
1360  else
1361  {
1362  nPerIndex.insert(collapseIndex, 1);
1363  }
1364  }
1365  }
1366 
1367  // 2. Size
1368  collapseStrings.resize(2*nPerIndex.size());
1369  forAllConstIter(Map<label>, nPerIndex, iter)
1370  {
1371  collapseStrings.insert(iter.key(), DynamicList<label>(iter()));
1372  }
1373 
1374  // 3. Fill
1375  forAll(allPointInfo, pointI)
1376  {
1377  const label collapseIndex = allPointInfo[pointI].collapseIndex();
1378 
1379  if (collapseIndex != -1 && collapseIndex != -2)
1380  {
1381  collapseStrings[collapseIndex].append(pointI);
1382  }
1383  }
1384  }
1385 
1386 
1387 
1388 
1389 // OFstream str2("collapseStrings_" + name(count) + ".obj");
1390 // // Dump point collapses
1391 // forAllConstIter(Map<DynamicList<label> >, collapseStrings, iter)
1392 // {
1393 // const label masterPoint = iter.key();
1394 // const DynamicList<label>& edgeCollapses = iter();
1395 //
1396 // forAll(edgeCollapses, eI)
1397 // {
1398 // meshTools::writeOBJ
1399 // (
1400 // str2,
1401 // mesh_.points()[edgeCollapses[eI]],
1402 // mesh_.points()[masterPoint]
1403 // );
1404 // }
1405 // }
1406 
1407 
1408 
1409  // Current faces (is also collapseStatus: f.size() < 3)
1410  faceList newFaces(mesh_.faces());
1411 
1412  // Current cellCollapse status
1413  boolList cellRemoved(mesh_.nCells(), false);
1414 
1415  label nUnvisited = 0;
1416  label nUncollapsed = 0;
1417  label nCollapsed = 0;
1418 
1419  forAll(allPointInfo, pI)
1420  {
1421  const pointEdgeCollapse& pec = allPointInfo[pI];
1422 
1423  if (pec.collapseIndex() == -1)
1424  {
1425  nUnvisited++;
1426  }
1427  else if (pec.collapseIndex() == -2)
1428  {
1429  nUncollapsed++;
1430  }
1431  else
1432  {
1433  nCollapsed++;
1434  }
1435  }
1436 
1437  label nPoints = allPointInfo.size();
1438 
1439  reduce(nPoints, sumOp<label>());
1440  reduce(nUnvisited, sumOp<label>());
1441  reduce(nUncollapsed, sumOp<label>());
1442  reduce(nCollapsed, sumOp<label>());
1443 
1444  Info<< incrIndent;
1445  Info<< indent << "Number of points : " << nPoints << nl
1446  << indent << "Not visited : " << nUnvisited << nl
1447  << indent << "Not collapsed : " << nUncollapsed << nl
1448  << indent << "Collapsed : " << nCollapsed << nl
1449  << endl;
1450  Info<< decrIndent;
1451 
1452  do
1453  {
1454  forAll(newFaces, faceI)
1455  {
1456  filterFace(collapseStrings, allPointInfo, newFaces[faceI]);
1457  }
1458 
1459  // Check if faces to be collapsed cause cells to become collapsed.
1460  label nCellCollapsed = 0;
1461 
1462  forAll(cells, cellI)
1463  {
1464  if (!cellRemoved[cellI])
1465  {
1466  const cell& cFaces = cells[cellI];
1467 
1468  label nFaces = cFaces.size();
1469 
1470  forAll(cFaces, i)
1471  {
1472  label faceI = cFaces[i];
1473 
1474  if (newFaces[faceI].size() < 3)
1475  {
1476  --nFaces;
1477 
1478  if (nFaces < 4)
1479  {
1480  Pout<< "Cell:" << cellI
1481  << " uses faces:" << cFaces
1482  << " of which too many are marked for removal:"
1483  << endl
1484  << " ";
1485 
1486 
1487  forAll(cFaces, j)
1488  {
1489  if (newFaces[cFaces[j]].size() < 3)
1490  {
1491  Pout<< ' '<< cFaces[j];
1492  }
1493  }
1494  Pout<< endl;
1495 
1496  cellRemoved[cellI] = true;
1497 
1498  // Collapse all edges of cell to nothing
1499 // collapseEdges(cellEdges[cellI]);
1500 
1501  nCellCollapsed++;
1502 
1503  break;
1504  }
1505  }
1506  }
1507  }
1508  }
1509 
1510  reduce(nCellCollapsed, sumOp<label>());
1511  Info<< indent << "Collapsing " << nCellCollapsed << " cells" << endl;
1512 
1513  if (nCellCollapsed == 0)
1514  {
1515  break;
1516  }
1517 
1518  } while (true);
1519 
1520 
1521  // Keep track of faces that have been done already.
1522  boolList doneFace(mesh_.nFaces(), false);
1523 
1524  {
1525  // Mark points used.
1526  boolList usedPoint(mesh_.nPoints(), false);
1527 
1528  forAll(cellRemoved, cellI)
1529  {
1530  if (cellRemoved[cellI])
1531  {
1532  meshMod.removeCell(cellI, -1);
1533  }
1534  }
1535 
1536  // Remove faces
1537  forAll(newFaces, faceI)
1538  {
1539  const face& f = newFaces[faceI];
1540 
1541  if (f.size() < 3)
1542  {
1543  meshMod.removeFace(faceI, -1);
1544  meshChanged = true;
1545 
1546  // Mark face as been done.
1547  doneFace[faceI] = true;
1548  }
1549  else
1550  {
1551  // Kept face. Mark vertices
1552  forAll(f, fp)
1553  {
1554  usedPoint[f[fp]] = true;
1555  }
1556  }
1557  }
1558 
1559  // Remove unused vertices that have not been marked for removal already
1560  forAll(usedPoint, pointI)
1561  {
1562  if (!usedPoint[pointI])
1563  {
1564  removedPoints[pointI] = true;
1565  meshMod.removePoint(pointI, -1);
1566  meshChanged = true;
1567  }
1568  }
1569  }
1570 
1571  // Modify the point location of the remaining points
1572  forAll(allPointInfo, pointI)
1573  {
1574  const label collapseIndex = allPointInfo[pointI].collapseIndex();
1575  const point& collapsePoint = allPointInfo[pointI].collapsePoint();
1576 
1577  if
1578  (
1579  removedPoints[pointI] == false
1580  && collapseIndex != -1
1581  && collapseIndex != -2
1582  )
1583  {
1584  meshMod.modifyPoint
1585  (
1586  pointI,
1587  collapsePoint,
1588  pointZones.whichZone(pointI),
1589  false
1590  );
1591  }
1592  }
1593 
1594 
1595  const polyBoundaryMesh& boundaryMesh = mesh_.boundaryMesh();
1596  const faceZoneMesh& faceZones = mesh_.faceZones();
1597 
1598  // Renumber faces that use points
1599  forAll(allPointInfo, pointI)
1600  {
1601  if (removedPoints[pointI] == true)
1602  {
1603  const labelList& changedFaces = pointFaces[pointI];
1604 
1605  forAll(changedFaces, changedFaceI)
1606  {
1607  label faceI = changedFaces[changedFaceI];
1608 
1609  if (!doneFace[faceI])
1610  {
1611  doneFace[faceI] = true;
1612 
1613  // Get current zone info
1614  label zoneID = faceZones.whichZone(faceI);
1615 
1616  bool zoneFlip = false;
1617 
1618  if (zoneID >= 0)
1619  {
1620  const faceZone& fZone = faceZones[zoneID];
1621 
1622  zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
1623  }
1624 
1625  // Get current connectivity
1626  label own = faceOwner[faceI];
1627  label nei = -1;
1628  label patchID = -1;
1629 
1630  if (mesh_.isInternalFace(faceI))
1631  {
1632  nei = faceNeighbour[faceI];
1633  }
1634  else
1635  {
1636  patchID = boundaryMesh.whichPatch(faceI);
1637  }
1638 
1639  meshMod.modifyFace
1640  (
1641  newFaces[faceI], // face
1642  faceI, // faceI to change
1643  own, // owner
1644  nei, // neighbour
1645  false, // flipFaceFlux
1646  patchID, // patch
1647  zoneID,
1648  zoneFlip
1649  );
1650 
1651  meshChanged = true;
1652  }
1653  }
1654  }
1655  }
1656 
1657  return meshChanged;
1658 }
1659 
1660 
1663  const globalIndex& globalPoints,
1664  const labelList& pointPriority,
1665  const Map<point>& collapsePointToLocation,
1666  PackedBoolList& collapseEdge,
1667  List<pointEdgeCollapse>& allPointInfo,
1668  const bool allowCellCollapse
1669 ) const
1670 {
1671  // Make sure we don't collapse cells
1672  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1673  const faceList faces = mesh_.faces();
1674  const edgeList& edges = mesh_.edges();
1675  const labelListList& faceEdges = mesh_.faceEdges();
1676  const labelListList& pointEdges = mesh_.pointEdges();
1677  const cellList& cells = mesh_.cells();
1678 
1679  labelHashSet uniqueCollapses;
1680  labelHashSet duplicateCollapses;
1681 
1682  while (true)
1683  {
1684  label nUncollapsed = 0;
1685 
1687  (
1688  mesh_,
1689  collapseEdge,
1691  0
1692  );
1693 
1694  // Create consistent set of collapses
1695  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1696  // Note: requires collapseEdge to be synchronised.
1697  syncCollapse
1698  (
1699  globalPoints,
1700  pointPriority,
1701  collapseEdge,
1702  collapsePointToLocation,
1703  allPointInfo
1704  );
1705 
1706  // Get collapsed faces
1707 
1708  PackedBoolList isCollapsedFace(mesh_.nFaces());
1709  PackedBoolList markedPoints(mesh_.nPoints());
1710 
1711  forAll(faces, faceI)
1712  {
1713  const face& f = faces[faceI];
1714 
1715  isCollapsedFace[faceI] = isFaceCollapsed(f, allPointInfo);
1716 
1717  if (isCollapsedFace[faceI] < 1)
1718  {
1719  determineDuplicatePointsOnFace
1720  (
1721  f,
1722  markedPoints,
1723  uniqueCollapses,
1724  duplicateCollapses,
1725  allPointInfo
1726  );
1727  }
1728  }
1729 
1730  // Synchronise the marked points
1732  (
1733  mesh_,
1734  markedPoints,
1736  0
1737  );
1738 
1739  // Mark all edges attached to the point for collapse
1740  forAll(markedPoints, pointI)
1741  {
1742  if (markedPoints[pointI])
1743  {
1744  const label index = allPointInfo[pointI].collapseIndex();
1745 
1746  const labelList& ptEdges = pointEdges[pointI];
1747 
1748  forAll(ptEdges, ptEdgeI)
1749  {
1750  const label edgeI = ptEdges[ptEdgeI];
1751  const label nbrPointI = edges[edgeI].otherVertex(pointI);
1752  const label nbrIndex
1753  = allPointInfo[nbrPointI].collapseIndex();
1754 
1755  if (collapseEdge[edgeI] && nbrIndex == index)
1756  {
1757  collapseEdge[edgeI] = false;
1758  nUncollapsed++;
1759  }
1760  }
1761  }
1762  }
1763 
1764  PackedBoolList markedEdges(mesh_.nEdges());
1765 
1766  if (!allowCellCollapse)
1767  {
1768  // Check collapsed cells
1769  forAll(cells, cellI)
1770  {
1771  const cell& cFaces = cells[cellI];
1772 
1773  label nFaces = cFaces.size();
1774 
1775  forAll(cFaces, fI)
1776  {
1777  label faceI = cFaces[fI];
1778 
1779  if (isCollapsedFace[faceI])
1780  {
1781  nFaces--;
1782  }
1783  }
1784 
1785  if (nFaces < 4)
1786  {
1787  forAll(cFaces, fI)
1788  {
1789  label faceI = cFaces[fI];
1790 
1791  const labelList& fEdges = faceEdges[faceI];
1792 
1793  // Unmark this face for collapse
1794  forAll(fEdges, fEdgeI)
1795  {
1796  label edgeI = fEdges[fEdgeI];
1797 
1798  if (collapseEdge[edgeI])
1799  {
1800  collapseEdge[edgeI] = false;
1801  nUncollapsed++;
1802  }
1803 
1804  markedEdges[edgeI] = true;
1805  }
1806 
1807  // Uncollapsed this face.
1808  isCollapsedFace[faceI] = false;
1809  nFaces++;
1810  }
1811  }
1812 
1813  if (nFaces < 4)
1814  {
1815  FatalErrorIn
1816  (
1817  "consistentCollapse\n"
1818  "(\n"
1819  " labelList&,\n"
1820  " Map<point>&,\n"
1821  " bool&,\n"
1822  ")"
1823  ) << "Cell " << cellI << " " << cFaces << nl
1824  << "is " << nFaces << ", "
1825  << "but cell collapse has been disabled."
1826  << abort(FatalError);
1827  }
1828  }
1829  }
1830 
1832  (
1833  mesh_,
1834  markedEdges,
1836  0
1837  );
1838 
1839  nUncollapsed += breakStringsAtEdges
1840  (
1841  markedEdges,
1842  collapseEdge,
1843  allPointInfo
1844  );
1845 
1846  reduce(nUncollapsed, sumOp<label>());
1847 
1848  Info<< " Uncollapsed edges = " << nUncollapsed << " / "
1849  << returnReduce(mesh_.nEdges(), sumOp<label>()) << endl;
1850 
1851  if (nUncollapsed == 0)
1852  {
1853  break;
1854  }
1855  }
1856 }
1857 
1858 
1861  const scalarField& minEdgeLen,
1862  const labelList& pointPriority,
1863  PackedBoolList& collapseEdge,
1864  Map<point>& collapsePointToLocation
1865 ) const
1866 {
1867  // Work out which edges to collapse
1868  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1869 
1870  const pointField& points = mesh_.points();
1871  const edgeList& edges = mesh_.edges();
1872 
1873  label nCollapsed = 0;
1874 
1875  forAll(edges, edgeI)
1876  {
1877  const edge& e = edges[edgeI];
1878 
1879  if (!collapseEdge[edgeI])
1880  {
1881  if (e.mag(points) < minEdgeLen[edgeI])
1882  {
1883  collapseEdge[edgeI] = true;
1884 
1885  label masterPointI = edgeMaster(pointPriority, e);
1886 
1887  if (masterPointI == -1)
1888  {
1889  const point average
1890  = 0.5*(points[e.start()] + points[e.end()]);
1891 
1892  collapsePointToLocation.set(e.start(), average);
1893  }
1894  else
1895  {
1896  const point& collapsePt = points[masterPointI];
1897 
1898  collapsePointToLocation.set(masterPointI, collapsePt);
1899  }
1900 
1901 
1902  nCollapsed++;
1903  }
1904  }
1905  }
1906 
1907  return nCollapsed;
1908 }
1909 
1910 
1913  const scalar maxCos,
1914  const labelList& pointPriority,
1915  PackedBoolList& collapseEdge,
1916  Map<point>& collapsePointToLocation
1917 ) const
1918 {
1919  const edgeList& edges = mesh_.edges();
1920  const pointField& points = mesh_.points();
1921  const labelListList& pointEdges = mesh_.pointEdges();
1922 
1923  // Point removal engine
1924  removePoints pointRemover(mesh_, false);
1925 
1926  // Find out points that can be deleted
1927  boolList pointCanBeDeleted;
1928  label nTotRemove = pointRemover.countPointUsage(maxCos, pointCanBeDeleted);
1929 
1930  // Rework point-to-remove into edge-to-collapse.
1931 
1932  label nCollapsed = 0;
1933 
1934  if (nTotRemove > 0)
1935  {
1936  forAll(pointEdges, pointI)
1937  {
1938  if (pointCanBeDeleted[pointI])
1939  {
1940  const labelList& pEdges = pointEdges[pointI];
1941 
1942  if (pEdges.size() == 2)
1943  {
1944  // Always the case?
1945 
1946  label e0 = pEdges[0];
1947  label e1 = pEdges[1];
1948 
1949  if (!collapseEdge[e0] && !collapseEdge[e1])
1950  {
1951  // Get lengths of both edges and choose the smallest
1952  scalar e0length = mag
1953  (
1954  points[edges[e0][0]] - points[edges[e0][1]]
1955  );
1956 
1957  scalar e1length = mag
1958  (
1959  points[edges[e1][0]] - points[edges[e1][1]]
1960  );
1961 
1962  if (e0length <= e1length)
1963  {
1964  collapseEdge[e0] = true;
1965 
1966  checkBoundaryPointMergeEdges
1967  (
1968  pointI,
1969  edges[e0].otherVertex(pointI),
1970  pointPriority,
1971  collapsePointToLocation
1972  );
1973  }
1974  else
1975  {
1976  collapseEdge[e1] = true;
1977 
1978  checkBoundaryPointMergeEdges
1979  (
1980  pointI,
1981  edges[e1].otherVertex(pointI),
1982  pointPriority,
1983  collapsePointToLocation
1984  );
1985  }
1986 
1987  nCollapsed++;
1988  }
1989  }
1990  }
1991  }
1992  }
1993 
1994  return nCollapsed;
1995 }
1996 
1997 
2000  const scalarField& faceFilterFactor,
2001  const labelList& pointPriority,
2002  PackedBoolList& collapseEdge,
2003  Map<point>& collapsePointToLocation
2004 ) const
2005 {
2006  const faceList& faces = mesh_.faces();
2007 
2008  const scalarField targetFaceSizes = calcTargetFaceSizes();
2009 
2010  // Calculate number of faces that will be collapsed to a point or an edge
2011  label nCollapseToPoint = 0;
2012  label nCollapseToEdge = 0;
2013 
2014  forAll(faces, fI)
2015  {
2016  const face& f = faces[fI];
2017 
2018  if (faceFilterFactor[fI] <= 0)
2019  {
2020  continue;
2021  }
2022 
2023  collapseType flagCollapseFace = collapseFace
2024  (
2025  pointPriority,
2026  f,
2027  fI,
2028  targetFaceSizes[fI],
2029  collapseEdge,
2030  collapsePointToLocation,
2031  faceFilterFactor
2032  );
2033 
2034  if (flagCollapseFace == noCollapse)
2035  {
2036  continue;
2037  }
2038  else if (flagCollapseFace == toPoint)
2039  {
2040  nCollapseToPoint++;
2041  }
2042  else if (flagCollapseFace == toEdge)
2043  {
2044  nCollapseToEdge++;
2045  }
2046  else
2047  {
2048  FatalErrorIn("collapseFaces(const polyMesh&, List<labelPair>&)")
2049  << "Face is marked to be collapsed to " << flagCollapseFace
2050  << ". Currently can only collapse to point/edge."
2051  << abort(FatalError);
2052  }
2053  }
2054 
2055  return labelPair(nCollapseToPoint, nCollapseToEdge);
2056 }
2057 
2058 
2061  const faceZone& fZone,
2062  const scalarField& faceFilterFactor,
2063  const labelList& pointPriority,
2064  PackedBoolList& collapseEdge,
2065  Map<point>& collapsePointToLocation
2066 ) const
2067 {
2068  const faceList& faces = mesh_.faces();
2069 
2070  const scalarField targetFaceSizes = calcTargetFaceSizes();
2071 
2072  // Calculate number of faces that will be collapsed to a point or an edge
2073  label nCollapseToPoint = 0;
2074  label nCollapseToEdge = 0;
2075 
2076  forAll(faces, fI)
2077  {
2078  if (fZone.whichFace(fI) == -1)
2079  {
2080  continue;
2081  }
2082 
2083  const face& f = faces[fI];
2084 
2085  if (faceFilterFactor[fI] <= 0)
2086  {
2087  continue;
2088  }
2089 
2090  collapseType flagCollapseFace = collapseFace
2091  (
2092  pointPriority,
2093  f,
2094  fI,
2095  targetFaceSizes[fI],
2096  collapseEdge,
2097  collapsePointToLocation,
2098  faceFilterFactor
2099  );
2100 
2101  if (flagCollapseFace == noCollapse)
2102  {
2103  continue;
2104  }
2105  else if (flagCollapseFace == toPoint)
2106  {
2107  nCollapseToPoint++;
2108  }
2109  else if (flagCollapseFace == toEdge)
2110  {
2111  nCollapseToEdge++;
2112  }
2113  else
2114  {
2115  FatalErrorIn("collapseFaces(const polyMesh&, List<labelPair>&)")
2116  << "Face is marked to be collapsed to " << flagCollapseFace
2117  << ". Currently can only collapse to point/edge."
2118  << abort(FatalError);
2119  }
2120  }
2121 
2122  return labelPair(nCollapseToPoint, nCollapseToEdge);
2123 
2124 // const edgeList& edges = mesh_.edges();
2125 // const pointField& points = mesh_.points();
2126 // const labelListList& edgeFaces = mesh_.edgeFaces();
2127 // const polyBoundaryMesh& bMesh = mesh_.boundaryMesh();
2128 //
2129 // forAll(edges, eI)
2130 // {
2131 // const edge& e = edges[eI];
2132 //
2133 // const labelList& eFaces = edgeFaces[eI];
2134 //
2135 // bool keepEdge = false;
2136 //
2137 // label nInternalFaces = 0;
2138 // label nPatchFaces = 0;
2139 // label nIndirectFaces = 0;
2140 //
2141 // bool coupled = false;
2142 //
2143 // forAll(eFaces, eFaceI)
2144 // {
2145 // const label eFaceIndex = eFaces[eFaceI];
2146 //
2147 // if (mesh_.isInternalFace(eFaceIndex))
2148 // {
2149 // nInternalFaces++;
2150 // }
2151 // else
2152 // {
2153 // const label patchIndex = bMesh.whichPatch(eFaceIndex);
2154 // const polyPatch& pPatch = bMesh[patchIndex];
2155 //
2156 // if (pPatch.coupled())
2157 // {
2158 // coupled = true;
2159 // nInternalFaces++;
2160 // }
2161 // else
2162 // {
2163 // // Keep the edge if an attached face is not in the zone
2164 // if (fZone.whichFace(eFaceIndex) == -1)
2165 // {
2166 // nPatchFaces++;
2167 // }
2168 // else
2169 // {
2170 // nIndirectFaces++;
2171 // }
2172 // }
2173 // }
2174 // }
2175 //
2176 // if (eFaces.size() != nInternalFaces + nPatchFaces + nIndirectFaces)
2177 // {
2178 // Pout<< eFaces.size() << " ("
2179 // << nInternalFaces << "/" << nPatchFaces << "/"
2180 // << nIndirectFaces << ")" << endl;
2181 // }
2182 //
2183 // if
2184 // (
2185 // eFaces.size() == nInternalFaces
2186 // || nIndirectFaces < (coupled ? 1 : 2)
2187 // )
2188 // {
2189 // keepEdge = true;
2190 // }
2191 //
2192 // if (!keepEdge)
2193 // {
2194 // collapseEdge[eI] = true;
2195 //
2196 // const Foam::point collapsePoint =
2197 // 0.5*(points[e.end()] + points[e.start()]);
2198 //
2199 // collapsePointToLocation.insert(e.start(), collapsePoint);
2200 // collapsePointToLocation.insert(e.end(), collapsePoint);
2201 // }
2202 // }
2203 
2204 // OFstream str
2205 // (
2206 // mesh_.time().path()
2207 // /"markedEdges_" + name(collapseEdge.count()) + ".obj"
2208 // );
2209 // label count = 0;
2210 //
2211 // forAll(collapseEdge, eI)
2212 // {
2213 // if (collapseEdge[eI])
2214 // {
2215 // const edge& e = edges[eI];
2216 //
2217 // meshTools::writeOBJ
2218 // (
2219 // str,
2220 // points[e.start()],
2221 // points[e.end()],
2222 // count
2223 // );
2224 // }
2225 // }
2226 }
2227 
2228 
2229 // ************************************************************************* //
const vectorField & faceAreas() const
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
label markMergeEdges(const scalar maxCos, const labelList &pointPriority, PackedBoolList &collapseEdge, Map< point > &collapsePointToLocation) const
Mark (in collapseEdge) any edges to merge.
label size() const
Return the number of elements in the VectorSpace = nCmpt.
Definition: VectorSpaceI.H:67
const pointField & points
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
bool setRefinement(const List< pointEdgeCollapse > &allPointInfo, polyTopoChange &meshMod) const
Play commands into polyTopoChange to create mesh.
label nFaces() const
dimensionedVector eigenValues(const dimensionedTensor &dt)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:310
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:232
dimensioned< scalar > mag(const dimensioned< Type > &)
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
const labelListList & pointEdges() const
labelList f(nPoints)
static label checkMeshQuality(const polyMesh &mesh, const dictionary &meshQualityDict, PackedBoolList &isErrorPoint)
Check mesh and mark points on faces in error.
Definition: edgeCollapser.C:82
PointFrompoint toPoint(const Foam::point &p)
bool empty() const
Return true if the UList is empty (ie, size() is zero).
Definition: UListI.H:313
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces)
Check mesh with mesh settings in dict. Collects incorrect faces.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
bool insert(const label &, const T &newElmt)
Insert a new hashedEntry.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
T & last()
Return the last element of the list.
Definition: UListI.H:131
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
static HashSet< label > checkBadFaces(const polyMesh &mesh, const dictionary &meshQualityDict)
Calls motionSmoother::checkMesh and returns a set of bad faces.
Definition: edgeCollapser.C:46
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:316
scalar mag(const pointField &) const
Magnitude of face area.
Definition: faceI.H:97
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.
A bit-packed bool list.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:242
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
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:310
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
label nCells() const
const cellList & cells() const
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:457
edgeList edges() const
Return edges in face point ordering,.
Definition: face.C:766
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh edges.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
label nEdges() const
messageStream Info
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void modifyPoint(const label, const point &, const label newZoneID, const bool inCell)
Modify coordinate.
Removes selected points from mesh and updates faces using these points.
Definition: removePoints.H:58
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
T & first()
Return the first element of the list.
Definition: UListI.H:117
point centre(const pointField &) const
Centre point of face.
Definition: face.C:493
const labelListList & faceEdges() const
Namespace for OpenFOAM.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:317
patches[0]
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: edgeI.H:163
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:379
static const char nl
Definition: Ostream.H:260
void setSize(const label)
Reset size of List.
Definition: List.C:318
const Cmpt & y() const
Definition: VectorI.H:71
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
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.
#define WarningIn(functionName)
Report a warning using Foam::Warning.
void clear()
Clear all entries from table.
Definition: HashTable.C:473
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1035
Foam::polyBoundaryMesh.
tensor inertia(const pointField &, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: face.C:731
dimensionedScalar det(const dimensionedSphericalTensor &dt)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
void removeCell(const label, const label)
Remove/merge cell.
Determines length of string of edges walked to point.
label otherVertex(const label a) const
Given one vertex, return the other.
Definition: edgeI.H:103
const Cmpt & x() const
Definition: VectorI.H:65
label markSmallEdges(const scalarField &minEdgeLen, const labelList &pointPriority, PackedBoolList &collapseEdge, Map< point > &collapsePointToLocation) const
Mark (in collapseEdge) any edges to collapse.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:429
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
static const label labelMin
Definition: label.H:61
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool set(const label &, const T &newElmt)
Assign a new hashedEntry, overwriting existing entries.
const labelListList & pointFaces() const
errorManip< error > abort(error &err)
Definition: errorManip.H:131
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1073
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
label longestEdge(const face &f, const pointField &pts)
Find the longest edge on a face. Face point labels index into pts.
Definition: face.C:875
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
const cellShapeList & cells
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
void modifyFace(const face &f, const label faceI, const label own, const label nei, const bool flipFaceFlux, const label patchID, const label zoneID, const bool zoneFlip)
Modify vertices or cell of face.
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
void removeFace(const label, const label)
Remove/merge face.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
label end() const
Return end vertex label.
Definition: edgeI.H:92
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
static const Vector zero
Definition: Vector.H:80
A HashTable with keys but without contents.
Definition: HashSet.H:59
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:85
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
Direct mesh changes based on v1.3 polyTopoChange syntax.
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:106
void removePoint(const label, const label)
Remove/merge point.
A List obtained as a section of another List.
Definition: SubList.H:53
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:97
label start() const
Return start vertex label.
Definition: edgeI.H:81
vector eigenVector(const tensor &, const scalar lambda)
Definition: tensor.C:207
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1060
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.
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:139
dimensionedSymmTensor sqr(const dimensionedVector &dv)
label nPoints
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1079
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:124
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
label nPoints() const
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
static const label labelMax
Definition: label.H:62
label countPointUsage(const scalar minCos, boolList &pointCanBeDeleted) const
Mark in pointCanBeDeleted the points that can be deleted.
Definition: removePoints.C:166
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116