primitiveMeshCheck.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-2023 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 "primitiveMesh.H"
27 #include "pyramidPointFaceRef.H"
28 #include "ListOps.H"
29 #include "unitConversion.H"
30 #include "SortableList.H"
31 #include "EdgeMap.H"
32 #include "primitiveMeshTools.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
37 (
38  const bool report
39 ) const
40 {
41  if (debug)
42  {
44  << "Checking whether the boundary is closed" << endl;
45  }
46 
47  const vectorField& areas = faceAreas();
48 
49  // Loop through all boundary faces and sum up the face area vectors.
50  // For a closed boundary, this should be zero in all vector components
51 
52  vector sumClosed(Zero);
53  scalar sumMagClosedBoundary = 0;
54 
55  for (label facei = nInternalFaces(); facei < areas.size(); facei++)
56  {
57  sumClosed += areas[facei];
58  sumMagClosedBoundary += mag(areas[facei]);
59  }
60 
61  reduce(sumClosed, sumOp<vector>());
62  reduce(sumMagClosedBoundary, sumOp<scalar>());
63 
64  vector openness = sumClosed/(sumMagClosedBoundary + vSmall);
65 
67  {
68  if (debug || report)
69  {
70  Info<< " ***Boundary openness " << openness
71  << " possible hole in boundary description."
72  << endl;
73  }
74 
75  return true;
76  }
77  else
78  {
79  if (debug || report)
80  {
81  Info<< " Boundary openness " << openness << " OK."
82  << endl;
83  }
84 
85  return false;
86  }
87 }
88 
89 
91 (
92  const bool report,
93  labelHashSet* setPtr,
94  labelHashSet* aspectSetPtr,
95  const Vector<label>& meshD
96 ) const
97 {
98  if (debug)
99  {
101  << "Checking whether cells are closed" << endl;
102  }
103 
104  const vectorField& faceAreas = this->faceAreas();
105  const scalarField& cellVolumes = this->cellVolumes();
106 
107  // Check that all cells labels are valid
108  const cellList& c = cells();
109 
110  label nErrorClosed = 0;
111 
112  forAll(c, cI)
113  {
114  const cell& curCell = c[cI];
115 
116  if (min(curCell) < 0 || max(curCell) > nFaces())
117  {
118  if (setPtr)
119  {
120  setPtr->insert(cI);
121  }
122 
123  nErrorClosed++;
124  }
125  }
126 
127  if (nErrorClosed > 0)
128  {
129  if (debug || report)
130  {
131  Info<< " ***Cells with invalid face labels found, number of cells "
132  << nErrorClosed << endl;
133  }
134 
135  return true;
136  }
137 
138 
139  scalarField openness;
140  scalarField aspectRatio;
142  (
143  *this,
144  meshD,
145  faceAreas,
146  cellVolumes,
147  openness,
148  aspectRatio
149  );
150 
151  label nOpen = 0;
152  scalar maxOpennessCell = max(openness);
153  label nAspect = 0;
154  scalar maxAspectRatio = max(aspectRatio);
155 
156  // Check the sums
157  forAll(openness, celli)
158  {
159  if (openness[celli] > polyMeshCheck::closedThreshold)
160  {
161  if (setPtr)
162  {
163  setPtr->insert(celli);
164  }
165 
166  nOpen++;
167  }
168 
169  if (aspectRatio[celli] > polyMeshCheck::aspectThreshold)
170  {
171  if (aspectSetPtr)
172  {
173  aspectSetPtr->insert(celli);
174  }
175 
176  nAspect++;
177  }
178  }
179 
180  reduce(nOpen, sumOp<label>());
181  reduce(maxOpennessCell, maxOp<scalar>());
182 
183  reduce(nAspect, sumOp<label>());
184  reduce(maxAspectRatio, maxOp<scalar>());
185 
186 
187  if (nOpen > 0)
188  {
189  if (debug || report)
190  {
191  Info<< " ***Open cells found, max cell openness: "
192  << maxOpennessCell << ", number of open cells " << nOpen
193  << endl;
194  }
195 
196  return true;
197  }
198 
199  if (nAspect > 0)
200  {
201  if (debug || report)
202  {
203  Info<< " ***High aspect ratio cells found, Max aspect ratio: "
204  << maxAspectRatio
205  << ", number of cells " << nAspect
206  << endl;
207  }
208 
209  return true;
210  }
211 
212  if (debug || report)
213  {
214  Info<< " Max cell openness = " << maxOpennessCell << " OK." << nl
215  << " Max aspect ratio = " << maxAspectRatio << " OK."
216  << endl;
217  }
218 
219  return false;
220 }
221 
222 
224 (
225  const bool report,
226  labelHashSet* setPtr
227 ) const
228 {
229  if (debug)
230  {
231  InfoInFunction << "Checking face area magnitudes" << endl;
232  }
233 
234  const vectorField& faceAreas = this->faceAreas();
235  const scalarField magFaceAreas(mag(faceAreas));
236 
237  scalar minArea = great;
238  scalar maxArea = -great;
239 
240  forAll(magFaceAreas, facei)
241  {
242  if (magFaceAreas[facei] < vSmall)
243  {
244  if (setPtr)
245  {
246  setPtr->insert(facei);
247  }
248  }
249 
250  minArea = min(minArea, magFaceAreas[facei]);
251  maxArea = max(maxArea, magFaceAreas[facei]);
252  }
253 
254  reduce(minArea, minOp<scalar>());
255  reduce(maxArea, maxOp<scalar>());
256 
257  if (minArea < vSmall)
258  {
259  if (debug || report)
260  {
261  Info<< " ***Zero or negative face area detected. "
262  "Minimum area: " << minArea << endl;
263  }
264 
265  return true;
266  }
267  else
268  {
269  if (debug || report)
270  {
271  Info<< " Minimum face area = " << minArea
272  << ". Maximum face area = " << maxArea
273  << ". Face area magnitudes OK." << endl;
274  }
275 
276  return false;
277  }
278 }
279 
280 
282 (
283  const bool report,
284  labelHashSet* setPtr
285 ) const
286 {
287  if (debug)
288  {
289  InfoInFunction << "Checking cell volumes" << endl;
290  }
291 
292  const scalarField& vols = cellVolumes();
293  scalar minVolume = great;
294  scalar maxVolume = -great;
295 
296  label nNegVolCells = 0;
297 
298  forAll(vols, celli)
299  {
300  if (vols[celli] < vSmall)
301  {
302  if (setPtr)
303  {
304  setPtr->insert(celli);
305  }
306 
307  nNegVolCells++;
308  }
309 
310  minVolume = min(minVolume, vols[celli]);
311  maxVolume = max(maxVolume, vols[celli]);
312  }
313 
314  reduce(minVolume, minOp<scalar>());
315  reduce(maxVolume, maxOp<scalar>());
316  reduce(nNegVolCells, sumOp<label>());
317 
318  if (minVolume < vSmall)
319  {
320  if (debug || report)
321  {
322  Info<< " ***Zero or negative cell volume detected. "
323  << "Minimum negative volume: " << minVolume
324  << ", Number of negative volume cells: " << nNegVolCells
325  << endl;
326  }
327 
328  return true;
329  }
330  else
331  {
332  if (debug || report)
333  {
334  Info<< " Min volume = " << minVolume
335  << ". Max volume = " << maxVolume
336  << ". Total volume = " << gSum(vols)
337  << ". Cell volumes OK." << endl;
338  }
339 
340  return false;
341  }
342 }
343 
344 
346 (
347  const bool report,
348  const scalar minPyrVol,
349  labelHashSet* setPtr
350 ) const
351 {
352  if (debug)
353  {
354  InfoInFunction << "Checking face orientation" << endl;
355  }
356 
357  const pointField& points = this->points();
358  const vectorField& ctrs = cellCentres();
359 
360 
361  scalarField ownPyrVol;
362  scalarField neiPyrVol;
364  (
365  *this,
366  points,
367  ctrs,
368  ownPyrVol,
369  neiPyrVol
370  );
371 
372 
373  label nErrorPyrs = 0;
374 
375  forAll(ownPyrVol, facei)
376  {
377  if (ownPyrVol[facei] < minPyrVol)
378  {
379  if (setPtr)
380  {
381  setPtr->insert(facei);
382  }
383 
384  nErrorPyrs++;
385  }
386 
387  if (isInternalFace(facei))
388  {
389  if (neiPyrVol[facei] < minPyrVol)
390  {
391  if (setPtr)
392  {
393  setPtr->insert(facei);
394  }
395  nErrorPyrs++;
396  }
397  }
398  }
399 
400  reduce(nErrorPyrs, sumOp<label>());
401 
402  if (nErrorPyrs > 0)
403  {
404  if (debug || report)
405  {
406  Info<< " ***Error in face pyramids: "
407  << nErrorPyrs << " faces are incorrectly oriented."
408  << endl;
409  }
410 
411  return true;
412  }
413  else
414  {
415  if (debug || report)
416  {
417  Info<< " Face pyramids OK." << endl;
418  }
419 
420  return false;
421  }
422 }
423 
424 
426 (
427  const bool report,
428  const scalar maxDeg,
429  labelHashSet* setPtr
430 ) const
431 {
432  if (debug)
433  {
434  InfoInFunction << "Checking face angles" << endl;
435  }
436 
437  if (maxDeg < -small || maxDeg > 180+small)
438  {
440  << "maxDeg should be [0..180] but is now " << maxDeg
441  << exit(FatalError);
442  }
443 
444  const scalar maxSin = Foam::sin(degToRad(maxDeg));
445 
446  const pointField& points = this->points();
447  const vectorField& faceAreas = this->faceAreas();
448 
450  (
451  maxSin,
452  *this,
453  points,
454  faceAreas
455  );
456  const scalarField& faceAngles = tfaceAngles();
457 
458  scalar maxEdgeSin = max(faceAngles);
459 
460  label nConcave = 0;
461 
462  forAll(faceAngles, facei)
463  {
464  if (faceAngles[facei] > small)
465  {
466  nConcave++;
467 
468  if (setPtr)
469  {
470  setPtr->insert(facei);
471  }
472  }
473  }
474 
475  reduce(nConcave, sumOp<label>());
476  reduce(maxEdgeSin, maxOp<scalar>());
477 
478  if (nConcave > 0)
479  {
480  scalar maxConcaveDegr =
481  radToDeg(Foam::asin(Foam::min(1.0, maxEdgeSin)));
482 
483  if (debug || report)
484  {
485  Info<< " *There are " << nConcave
486  << " faces with concave angles between consecutive"
487  << " edges. Max concave angle = " << maxConcaveDegr
488  << " degrees." << endl;
489  }
490 
491  return true;
492  }
493  else
494  {
495  if (debug || report)
496  {
497  Info<< " All angles in faces OK." << endl;
498  }
499 
500  return false;
501  }
502 }
503 
504 
506 (
507  const bool report,
508  const scalar warnFlatness,
509  labelHashSet* setPtr
510 ) const
511 {
512  if (debug)
513  {
514  InfoInFunction << "Checking face flatness" << endl;
515  }
516 
517  if (warnFlatness < 0 || warnFlatness > 1)
518  {
520  << "warnFlatness should be [0..1] but is now " << warnFlatness
521  << exit(FatalError);
522  }
523 
524  const pointField& points = this->points();
525  const vectorField& faceCentres = this->faceCentres();
526  const vectorField& faceAreas = this->faceAreas();
527  const faceList& fcs = faces();
528 
530  (
531  *this,
532  points,
533  faceCentres,
534  faceAreas
535  );
536  const scalarField& faceFlatness = tfaceFlatness();
537 
538  scalarField magAreas(mag(faceAreas));
539 
540  scalar minFlatness = great;
541  scalar sumFlatness = 0;
542  label nSummed = 0;
543  label nWarped = 0;
544 
545  forAll(faceFlatness, facei)
546  {
547  if (fcs[facei].size() > 3 && magAreas[facei] > vSmall)
548  {
549  sumFlatness += faceFlatness[facei];
550  nSummed++;
551 
552  minFlatness = min(minFlatness, faceFlatness[facei]);
553 
554  if (faceFlatness[facei] < warnFlatness)
555  {
556  nWarped++;
557 
558  if (setPtr)
559  {
560  setPtr->insert(facei);
561  }
562  }
563  }
564  }
565 
566 
567  reduce(nWarped, sumOp<label>());
568  reduce(minFlatness, minOp<scalar>());
569 
570  reduce(nSummed, sumOp<label>());
571  reduce(sumFlatness, sumOp<scalar>());
572 
573  if (debug || report)
574  {
575  if (nSummed > 0)
576  {
577  Info<< " Face flatness (1 = flat, 0 = butterfly) : min = "
578  << minFlatness << " average = " << sumFlatness / nSummed
579  << endl;
580  }
581  }
582 
583 
584  if (nWarped> 0)
585  {
586  if (debug || report)
587  {
588  Info<< " *There are " << nWarped
589  << " faces with ratio between projected and actual area < "
590  << warnFlatness << endl;
591 
592  Info<< " Minimum ratio (minimum flatness, maximum warpage) = "
593  << minFlatness << endl;
594  }
595 
596  return true;
597  }
598  else
599  {
600  if (debug || report)
601  {
602  Info<< " All face flatness OK." << endl;
603  }
604 
605  return false;
606  }
607 }
608 
609 
611 (
612  const bool report,
613  labelHashSet* setPtr
614 ) const
615 {
616  if (debug)
617  {
618  InfoInFunction << "Checking for concave cells" << endl;
619  }
620 
621  const vectorField& fAreas = faceAreas();
622  const pointField& fCentres = faceCentres();
623  const cellList& c = cells();
624  const labelList& fOwner = faceOwner();
625 
626  label nConcaveCells = 0;
627 
628  forAll(c, celli)
629  {
630  const cell& cFaces = c[celli];
631 
632  bool concave = false;
633 
634  forAll(cFaces, i)
635  {
636  if (concave)
637  {
638  break;
639  }
640 
641  label fI = cFaces[i];
642 
643  const point& fC = fCentres[fI];
644 
645  vector fN = fAreas[fI];
646 
647  fN /= max(mag(fN), vSmall);
648 
649  // Flip normal if required so that it is always pointing out of
650  // the cell
651  if (fOwner[fI] != celli)
652  {
653  fN *= -1;
654  }
655 
656  // Is the centre of any other face of the cell on the
657  // wrong side of the plane of this face?
658 
659  forAll(cFaces, j)
660  {
661  if (j != i)
662  {
663  label fJ = cFaces[j];
664 
665  const point& pt = fCentres[fJ];
666 
667  // If the cell is concave, the point will be on the
668  // positive normal side of the plane of f, defined by
669  // its centre and normal, and the angle between (pt -
670  // fC) and fN will be less than 90 degrees, so the dot
671  // product will be positive.
672 
673  vector pC = (pt - fC);
674 
675  pC /= max(mag(pC), vSmall);
676 
677  if ((pC & fN) > -polyMeshCheck::planarCosAngle)
678  {
679  // Concave or planar face
680 
681  concave = true;
682 
683  if (setPtr)
684  {
685  setPtr->insert(celli);
686  }
687 
688  nConcaveCells++;
689 
690  break;
691  }
692  }
693  }
694  }
695  }
696 
697  reduce(nConcaveCells, sumOp<label>());
698 
699  if (nConcaveCells > 0)
700  {
701  if (debug || report)
702  {
703  Info<< " ***Concave cells (using face planes) found,"
704  << " number of cells: " << nConcaveCells << endl;
705  }
706 
707  return true;
708  }
709  else
710  {
711  if (debug || report)
712  {
713  Info<< " Concave cell check OK." << endl;
714  }
715 
716  return false;
717  }
718 
719  return false;
720 }
721 
722 
724 (
725  const bool report,
726  labelHashSet* setPtr
727 ) const
728 {
729  if (debug)
730  {
731  InfoInFunction << "Checking face ordering" << endl;
732  }
733 
734  // Check whether internal faces are ordered in the upper triangular order
735  const labelList& own = faceOwner();
736  const labelList& nei = faceNeighbour();
737 
738  const cellList& c = cells();
739 
740  label internal = nInternalFaces();
741 
742  // Has error occurred?
743  bool error = false;
744  // Have multiple faces been detected?
745  label nMultipleCells = false;
746 
747  // Loop through faceCells once more and make sure that for internal cell
748  // the first label is smaller
749  for (label facei = 0; facei < internal; facei++)
750  {
751  if (own[facei] >= nei[facei])
752  {
753  error = true;
754 
755  if (setPtr)
756  {
757  setPtr->insert(facei);
758  }
759  }
760  }
761 
762  // Loop through all cells. For each cell, find the face that is internal
763  // and add it to the check list (upper triangular order).
764  // Once the list is completed, check it against the faceCell list
765 
766  forAll(c, celli)
767  {
768  const labelList& curFaces = c[celli];
769 
770  // Neighbouring cells
771  SortableList<label> nbr(curFaces.size());
772 
773  forAll(curFaces, i)
774  {
775  label facei = curFaces[i];
776 
777  if (facei >= nInternalFaces())
778  {
779  // Sort last
780  nbr[i] = labelMax;
781  }
782  else
783  {
784  label nbrCelli = nei[facei];
785 
786  if (nbrCelli == celli)
787  {
788  nbrCelli = own[facei];
789  }
790 
791  if (celli < nbrCelli)
792  {
793  // celli is master
794  nbr[i] = nbrCelli;
795  }
796  else
797  {
798  // nbrCell is master. Let it handle this face.
799  nbr[i] = labelMax;
800  }
801  }
802  }
803 
804  nbr.sort();
805 
806  // Now nbr holds the cellCells in incremental order. Check:
807  // - neighbouring cells appear only once. Since nbr is sorted this
808  // is simple check on consecutive elements
809  // - faces indexed in same order as nbr are incrementing as well.
810 
811  label prevCell = nbr[0];
812  label prevFace = curFaces[nbr.indices()[0]];
813 
814  bool hasMultipleFaces = false;
815 
816  for (label i = 1; i < nbr.size(); i++)
817  {
818  label thisCell = nbr[i];
819  label thisFace = curFaces[nbr.indices()[i]];
820 
821  if (thisCell == labelMax)
822  {
823  break;
824  }
825 
826  if (thisCell == prevCell)
827  {
828  hasMultipleFaces = true;
829 
830  if (setPtr)
831  {
832  setPtr->insert(prevFace);
833  setPtr->insert(thisFace);
834  }
835  }
836  else if (thisFace < prevFace)
837  {
838  error = true;
839 
840  if (setPtr)
841  {
842  setPtr->insert(thisFace);
843  }
844  }
845 
846  prevCell = thisCell;
847  prevFace = thisFace;
848  }
849 
850  if (hasMultipleFaces)
851  {
852  nMultipleCells++;
853  }
854  }
855 
856  reduce(error, orOp<bool>());
857  reduce(nMultipleCells, sumOp<label>());
858 
859  if ((debug || report) && nMultipleCells > 0)
860  {
861  Info<< " <<Found " << nMultipleCells
862  << " neighbouring cells with multiple in between faces." << endl;
863  }
864 
865  if (error)
866  {
867  if (debug || report)
868  {
869  Info<< " ***Faces not in upper triangular order." << endl;
870  }
871 
872  return true;
873  }
874  else
875  {
876  if (debug || report)
877  {
878  Info<< " Upper triangular ordering OK." << endl;
879  }
880 
881  return false;
882  }
883 }
884 
885 
887 (
888  const bool report,
889  labelHashSet* setPtr
890 ) const
891 {
892  if (debug)
893  {
894  InfoInFunction << "Checking topological cell openness" << endl;
895  }
896 
897  label nOpenCells = 0;
898 
899  const faceList& f = faces();
900  const cellList& c = cells();
901 
902  forAll(c, celli)
903  {
904  const labelList& curFaces = c[celli];
905 
906  const edgeList cellEdges = c[celli].edges(f);
907 
908  labelList edgeUsage(cellEdges.size(), 0);
909 
910  forAll(curFaces, facei)
911  {
912  edgeList curFaceEdges = f[curFaces[facei]].edges();
913 
914  forAll(curFaceEdges, faceEdgeI)
915  {
916  const edge& curEdge = curFaceEdges[faceEdgeI];
917 
918  forAll(cellEdges, cellEdgeI)
919  {
920  if (cellEdges[cellEdgeI] == curEdge)
921  {
922  edgeUsage[cellEdgeI]++;
923  break;
924  }
925  }
926  }
927  }
928 
929  edgeList singleEdges(cellEdges.size());
930  label nSingleEdges = 0;
931 
932  forAll(edgeUsage, edgeI)
933  {
934  if (edgeUsage[edgeI] == 1)
935  {
936  singleEdges[nSingleEdges] = cellEdges[edgeI];
937  nSingleEdges++;
938  }
939  else if (edgeUsage[edgeI] != 2)
940  {
941  if (setPtr)
942  {
943  setPtr->insert(celli);
944  }
945  }
946  }
947 
948  if (nSingleEdges > 0)
949  {
950  if (setPtr)
951  {
952  setPtr->insert(celli);
953  }
954 
955  nOpenCells++;
956  }
957  }
958 
959  reduce(nOpenCells, sumOp<label>());
960 
961  if (nOpenCells > 0)
962  {
963  if (debug || report)
964  {
965  Info<< " ***Open cells found, number of cells: " << nOpenCells
966  << ". This problem may be fixable using the zipUpMesh utility."
967  << endl;
968  }
969 
970  return true;
971  }
972  else
973  {
974  if (debug || report)
975  {
976  Info<< " Topological cell zip-up check OK." << endl;
977  }
978 
979  return false;
980  }
981 }
982 
983 
985 (
986  const bool report,
987  labelHashSet* setPtr
988 ) const
989 {
990  if (debug)
991  {
992  InfoInFunction << "Checking face vertices" << endl;
993  }
994 
995  // Check that all vertex labels are valid
996  const faceList& f = faces();
997 
998  label nErrorFaces = 0;
999 
1000  forAll(f, fI)
1001  {
1002  const face& curFace = f[fI];
1003 
1004  if (min(curFace) < 0 || max(curFace) > nPoints())
1005  {
1006  if (setPtr)
1007  {
1008  setPtr->insert(fI);
1009  }
1010 
1011  nErrorFaces++;
1012  }
1013 
1014  // Uniqueness of vertices
1015  labelHashSet facePoints(2*curFace.size());
1016 
1017  forAll(curFace, fp)
1018  {
1019  bool inserted = facePoints.insert(curFace[fp]);
1020 
1021  if (!inserted)
1022  {
1023  if (setPtr)
1024  {
1025  setPtr->insert(fI);
1026  }
1027 
1028  nErrorFaces++;
1029  }
1030  }
1031  }
1032 
1033  reduce(nErrorFaces, sumOp<label>());
1034 
1035  if (nErrorFaces > 0)
1036  {
1037  if (debug || report)
1038  {
1039  Info<< " Faces with invalid vertex labels found, "
1040  << " number of faces: " << nErrorFaces << endl;
1041  }
1042 
1043  return true;
1044  }
1045  else
1046  {
1047  if (debug || report)
1048  {
1049  Info<< " Face vertices OK." << endl;
1050  }
1051 
1052  return false;
1053  }
1054 }
1055 
1056 
1058 (
1059  const bool report,
1060  labelHashSet* setPtr
1061 ) const
1062 {
1063  if (debug)
1064  {
1065  InfoInFunction << "Checking points" << endl;
1066  }
1067 
1068  label nFaceErrors = 0;
1069  label nCellErrors = 0;
1070 
1071  const labelListList& pf = pointFaces();
1072 
1073  forAll(pf, pointi)
1074  {
1075  if (pf[pointi].empty())
1076  {
1077  if (setPtr)
1078  {
1079  setPtr->insert(pointi);
1080  }
1081 
1082  nFaceErrors++;
1083  }
1084  }
1085 
1086 
1087  forAll(pf, pointi)
1088  {
1089  const labelList& pc = pointCells(pointi);
1090 
1091  if (pc.empty())
1092  {
1093  if (setPtr)
1094  {
1095  setPtr->insert(pointi);
1096  }
1097 
1098  nCellErrors++;
1099  }
1100  }
1101 
1102  reduce(nFaceErrors, sumOp<label>());
1103  reduce(nCellErrors, sumOp<label>());
1104 
1105  if (nFaceErrors > 0 || nCellErrors > 0)
1106  {
1107  if (debug || report)
1108  {
1109  Info<< " ***Unused points found in the mesh, "
1110  "number unused by faces: " << nFaceErrors
1111  << " number unused by cells: " << nCellErrors
1112  << endl;
1113  }
1114 
1115  return true;
1116  }
1117  else
1118  {
1119  if (debug || report)
1120  {
1121  Info<< " Point usage OK." << endl;
1122  }
1123 
1124  return false;
1125  }
1126 }
1127 
1128 
1130 (
1131  const label facei,
1132  const Map<label>& nCommonPoints,
1133  label& nBaffleFaces,
1134  labelHashSet* setPtr
1135 ) const
1136 {
1137  bool error = false;
1138 
1139  forAllConstIter(Map<label>, nCommonPoints, iter)
1140  {
1141  label nbFacei = iter.key();
1142  label nCommon = iter();
1143 
1144  const face& curFace = faces()[facei];
1145  const face& nbFace = faces()[nbFacei];
1146 
1147  if (nCommon == nbFace.size() || nCommon == curFace.size())
1148  {
1149  if (nbFace.size() != curFace.size())
1150  {
1151  error = true;
1152  }
1153  else
1154  {
1155  nBaffleFaces++;
1156  }
1157 
1158  if (setPtr)
1159  {
1160  setPtr->insert(facei);
1161  setPtr->insert(nbFacei);
1162  }
1163  }
1164  }
1165 
1166  return error;
1167 }
1168 
1169 
1171 (
1172  const label facei,
1173  const Map<label>& nCommonPoints,
1174  labelHashSet* setPtr
1175 ) const
1176 {
1177  bool error = false;
1178 
1179  forAllConstIter(Map<label>, nCommonPoints, iter)
1180  {
1181  label nbFacei = iter.key();
1182  label nCommon = iter();
1183 
1184  const face& curFace = faces()[facei];
1185  const face& nbFace = faces()[nbFacei];
1186 
1187  if
1188  (
1189  nCommon >= 2
1190  && nCommon != nbFace.size()
1191  && nCommon != curFace.size()
1192  )
1193  {
1194  forAll(curFace, fp)
1195  {
1196  // Get the index in the neighbouring face shared with curFace
1197  label nb = findIndex(nbFace, curFace[fp]);
1198 
1199  if (nb != -1)
1200  {
1201 
1202  // Check the whole face from nb onwards for shared vertices
1203  // with neighbouring face. Rule is that any shared vertices
1204  // should be consecutive on both faces i.e. if they are
1205  // vertices fp,fp+1,fp+2 on one face they should be
1206  // vertices nb, nb+1, nb+2 (or nb+2, nb+1, nb) on the
1207  // other face.
1208 
1209 
1210  // Vertices before and after on curFace
1211  label fpPlus1 = curFace.fcIndex(fp);
1212  label fpMin1 = curFace.rcIndex(fp);
1213 
1214  // Vertices before and after on nbFace
1215  label nbPlus1 = nbFace.fcIndex(nb);
1216  label nbMin1 = nbFace.rcIndex(nb);
1217 
1218  // Find order of walking by comparing next points on both
1219  // faces.
1220  label curInc = labelMax;
1221  label nbInc = labelMax;
1222 
1223  if (nbFace[nbPlus1] == curFace[fpPlus1])
1224  {
1225  curInc = 1;
1226  nbInc = 1;
1227  }
1228  else if (nbFace[nbPlus1] == curFace[fpMin1])
1229  {
1230  curInc = -1;
1231  nbInc = 1;
1232  }
1233  else if (nbFace[nbMin1] == curFace[fpMin1])
1234  {
1235  curInc = -1;
1236  nbInc = -1;
1237  }
1238  else
1239  {
1240  curInc = 1;
1241  nbInc = -1;
1242  }
1243 
1244 
1245  // Pass1: loop until start of common vertices found.
1246  label curNb = nb;
1247  label curFp = fp;
1248 
1249  do
1250  {
1251  curFp += curInc;
1252 
1253  if (curFp >= curFace.size())
1254  {
1255  curFp = 0;
1256  }
1257  else if (curFp < 0)
1258  {
1259  curFp = curFace.size()-1;
1260  }
1261 
1262  curNb += nbInc;
1263 
1264  if (curNb >= nbFace.size())
1265  {
1266  curNb = 0;
1267  }
1268  else if (curNb < 0)
1269  {
1270  curNb = nbFace.size()-1;
1271  }
1272  } while (curFace[curFp] == nbFace[curNb]);
1273 
1274 
1275  // Pass2: check equality walking from curFp, curNb
1276  // in opposite order.
1277 
1278  curInc = -curInc;
1279  nbInc = -nbInc;
1280 
1281  for (label commonI = 0; commonI < nCommon; commonI++)
1282  {
1283  curFp += curInc;
1284 
1285  if (curFp >= curFace.size())
1286  {
1287  curFp = 0;
1288  }
1289  else if (curFp < 0)
1290  {
1291  curFp = curFace.size()-1;
1292  }
1293 
1294  curNb += nbInc;
1295 
1296  if (curNb >= nbFace.size())
1297  {
1298  curNb = 0;
1299  }
1300  else if (curNb < 0)
1301  {
1302  curNb = nbFace.size()-1;
1303  }
1304 
1305  if (curFace[curFp] != nbFace[curNb])
1306  {
1307  if (setPtr)
1308  {
1309  setPtr->insert(facei);
1310  setPtr->insert(nbFacei);
1311  }
1312 
1313  error = true;
1314 
1315  break;
1316  }
1317  }
1318 
1319 
1320  // Done the curFace - nbFace combination.
1321  break;
1322  }
1323  }
1324  }
1325  }
1326 
1327  return error;
1328 }
1329 
1330 
1332 (
1333  const bool report,
1334  labelHashSet* setPtr
1335 ) const
1336 {
1337  if (debug)
1338  {
1339  InfoInFunction << "Checking face-face connectivity" << endl;
1340  }
1341 
1342  const labelListList& pf = pointFaces();
1343 
1344  label nBaffleFaces = 0;
1345  label nErrorDuplicate = 0;
1346  label nErrorOrder = 0;
1347  Map<label> nCommonPoints(100);
1348 
1349  for (label facei = 0; facei < nFaces(); facei++)
1350  {
1351  const face& curFace = faces()[facei];
1352 
1353  // Calculate number of common points between current facei and
1354  // neighbouring face. Store on map.
1355  nCommonPoints.clear();
1356 
1357  forAll(curFace, fp)
1358  {
1359  label pointi = curFace[fp];
1360 
1361  const labelList& nbs = pf[pointi];
1362 
1363  forAll(nbs, nbI)
1364  {
1365  label nbFacei = nbs[nbI];
1366 
1367  if (facei < nbFacei)
1368  {
1369  // Only check once for each combination of two faces.
1370 
1371  Map<label>::iterator fnd = nCommonPoints.find(nbFacei);
1372 
1373  if (fnd == nCommonPoints.end())
1374  {
1375  // First common vertex found.
1376  nCommonPoints.insert(nbFacei, 1);
1377  }
1378  else
1379  {
1380  fnd()++;
1381  }
1382  }
1383  }
1384  }
1385 
1386  // Perform various checks on common points
1387 
1388  // Check all vertices shared (duplicate point)
1389  if (checkDuplicateFaces(facei, nCommonPoints, nBaffleFaces, setPtr))
1390  {
1391  nErrorDuplicate++;
1392  }
1393 
1394  // Check common vertices are consecutive on both faces
1395  if (checkCommonOrder(facei, nCommonPoints, setPtr))
1396  {
1397  nErrorOrder++;
1398  }
1399  }
1400 
1401  reduce(nBaffleFaces, sumOp<label>());
1402  reduce(nErrorDuplicate, sumOp<label>());
1403  reduce(nErrorOrder, sumOp<label>());
1404 
1405  if (nBaffleFaces)
1406  {
1407  Info<< " Number of identical duplicate faces (baffle faces): "
1408  << nBaffleFaces << endl;
1409  }
1410 
1411  if (nErrorDuplicate > 0 || nErrorOrder > 0)
1412  {
1413  // These are actually warnings, not errors.
1414  if (nErrorDuplicate > 0)
1415  {
1416  Info<< " <<Number of duplicate (not baffle) faces found: "
1417  << nErrorDuplicate
1418  << ". This might indicate a problem." << endl;
1419  }
1420 
1421  if (nErrorOrder > 0)
1422  {
1423  Info<< " <<Number of faces with non-consecutive shared points: "
1424  << nErrorOrder << ". This might indicate a problem." << endl;
1425  }
1426 
1427  return false; // return true;
1428  }
1429  else
1430  {
1431  if (debug || report)
1432  {
1433  Info<< " Face-face connectivity OK." << endl;
1434  }
1435 
1436  return false;
1437  }
1438 }
1439 
1440 
1441 // ************************************************************************* //
Various functions to operate on Lists.
#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
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
void clear()
Clear all entries from table.
Definition: HashTable.C:468
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: SortableList.H:55
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:96
void sort()
(stable) sort the list (if changed after construction time)
Definition: SortableList.C:112
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:71
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
static tmp< scalarField > faceConcavity(const scalar maxSin, const primitiveMesh &mesh, const pointField &p, const vectorField &faceAreas)
Generate face concavity field. Returns per face the (sin of the)
static void cellClosedness(const primitiveMesh &mesh, const Vector< label > &meshD, const vectorField &areas, const scalarField &vols, scalarField &openness, scalarField &aratio)
Generate cell openness and cell aspect ratio field.
static void facePyramidVolume(const primitiveMesh &mesh, const pointField &points, const vectorField &cellCtrs, scalarField &ownPyrVol, scalarField &neiPyrVol)
Generate face pyramid volume fields.
static tmp< scalarField > faceFlatness(const primitiveMesh &mesh, const pointField &p, const vectorField &fCtrs, const vectorField &faceAreas)
Generate face flatness field. Compares the individual triangles'.
bool checkClosedCells(const bool report=false, labelHashSet *setPtr=nullptr, labelHashSet *highAspectSetPtr=nullptr, const Vector< label > &solutionD=Vector< label >::one) const
Check cells for closedness.
bool checkFacePyramids(const bool report=false, const scalar minPyrVol=-small, labelHashSet *setPtr=nullptr) const
Check face pyramid volume.
bool checkCellsZipUp(const bool report=false, labelHashSet *setPtr=nullptr) const
Check cell zip-up.
bool checkFaceAreas(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for negative face areas.
bool checkCommonOrder(const label, const Map< label > &, labelHashSet *) const
Check that shared points are in consecutive order.
bool checkFaceAngles(const bool report=false, const scalar maxSin=10, labelHashSet *setPtr=nullptr) const
Check face angles.
bool checkUpperTriangular(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face ordering.
bool checkConcaveCells(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for concave cells by the planes of faces.
bool checkClosedBoundary(const bool report=false) const
Check boundary for closedness.
bool checkCellVolumes(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for negative cell volumes.
bool checkDuplicateFaces(const label, const Map< label > &, label &nBaffleFaces, labelHashSet *) const
Check if all points on face are shared with another face.
bool checkFaceVertices(const bool report=false, labelHashSet *setPtr=nullptr) const
Check uniqueness of face vertices.
label nInternalFaces() const
bool checkPoints(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for unused points.
const vectorField & faceAreas() const
bool checkFaceFlatness(const bool report, const scalar warnFlatness, labelHashSet *setPtr) const
Check face warpage: decompose face and check ratio between.
bool checkFaceFaces(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face-face connectivity.
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const pointField & points
label nPoints
const cellShapeList & cells
#define InfoInFunction
Report an information message using Foam::Info.
const dimensionedScalar c
Speed of light in a vacuum.
scalar closedThreshold
Data to control mesh checking.
Definition: polyMeshCheck.C:33
scalar aspectThreshold
Aspect ratio warning threshold.
Definition: polyMeshCheck.C:34
scalar planarCosAngle
Threshold where faces are considered coplanar.
Definition: polyMeshCheck.C:37
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
dimensionedScalar asin(const dimensionedScalar &ds)
Type gSum(const FieldField< Field, Type > &f)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:275
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
static const label labelMax
Definition: label.H:62
static const char nl
Definition: Ostream.H:260
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
labelList f(nPoints)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
Unit conversion functions.