primitiveMeshCheck.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "primitiveMesh.H"
30 #include "pyramidPointFaceRef.H"
31 #include "ListOps.H"
32 #include "unitConversion.H"
33 #include "SortableList.H"
34 #include "edgeHashes.H"
35 #include "primitiveMeshTools.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 Foam::scalar Foam::primitiveMesh::closedThreshold_ = 1.0e-6;
40 Foam::scalar Foam::primitiveMesh::aspectThreshold_ = 1000;
41 Foam::scalar Foam::primitiveMesh::nonOrthThreshold_ = 70; // deg
43 Foam::scalar Foam::primitiveMesh::planarCosAngle_ = 1.0e-6;
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
49 (
50  const vectorField& areas,
51  const bool report,
52  const bitSet& internalOrCoupledFaces
53 ) const
54 {
55  DebugInFunction << "Checking if boundary is closed" << endl;
56 
57  // Loop through all boundary faces and sum up the face area vectors.
58  // For a closed boundary, this should be zero in all vector components
59 
60  Vector<solveScalar> sumClosed(Zero);
61  solveScalar sumMagClosedBoundary = 0;
62 
63  for (label facei = nInternalFaces(); facei < areas.size(); facei++)
64  {
65  if (!internalOrCoupledFaces.size() || !internalOrCoupledFaces[facei])
66  {
67  sumClosed += Vector<solveScalar>(areas[facei]);
68  sumMagClosedBoundary += mag(areas[facei]);
69  }
70  }
71 
72  reduce(sumClosed, sumOp<Vector<solveScalar>>());
73  reduce(sumMagClosedBoundary, sumOp<solveScalar>());
74 
75  Vector<solveScalar> openness = sumClosed/(sumMagClosedBoundary + VSMALL);
76 
77  if (cmptMax(cmptMag(openness)) > closedThreshold_)
78  {
79  if (debug || report)
80  {
81  Info<< " ***Boundary openness " << openness
82  << " possible hole in boundary description."
83  << endl;
84  }
85 
86  return true;
87  }
88 
89  if (debug || report)
90  {
91  Info<< " Boundary openness " << openness << " OK."
92  << endl;
93  }
94 
95  return false;
96 }
97 
98 
100 (
101  const vectorField& faceAreas,
102  const scalarField& cellVolumes,
103  const bool report,
104  labelHashSet* setPtr,
105  labelHashSet* aspectSetPtr,
106  const Vector<label>& meshD
107 ) const
108 {
109  DebugInFunction << "Checking if cells are closed" << endl;
110 
111  // Check that all cells labels are valid
112  const cellList& c = cells();
113 
114  label nErrorClosed = 0;
115 
116  forAll(c, cI)
117  {
118  const cell& curCell = c[cI];
119 
120  if (min(curCell) < 0 || max(curCell) > nFaces())
121  {
122  if (setPtr)
123  {
124  setPtr->insert(cI);
125  }
126 
127  nErrorClosed++;
128  }
129  }
130 
131  if (nErrorClosed > 0)
132  {
133  if (debug || report)
134  {
135  Info<< " ***Cells with invalid face labels found, number of cells "
136  << nErrorClosed << endl;
137  }
138 
139  return true;
140  }
141 
142 
143  scalarField openness;
144  scalarField aspectRatio;
145  primitiveMeshTools::cellClosedness
146  (
147  *this,
148  meshD,
149  faceAreas,
150  cellVolumes,
151  openness,
152  aspectRatio
153  );
154 
155  label nOpen = 0;
156  scalar maxOpennessCell = max(openness);
157  label nAspect = 0;
158  scalar maxAspectRatio = max(aspectRatio);
159 
160  // Check the sums
161  forAll(openness, celli)
162  {
163  if (openness[celli] > closedThreshold_)
164  {
165  if (setPtr)
166  {
167  setPtr->insert(celli);
168  }
169 
170  nOpen++;
171  }
172 
173  if (aspectRatio[celli] > aspectThreshold_)
174  {
175  if (aspectSetPtr)
176  {
177  aspectSetPtr->insert(celli);
178  }
179 
180  nAspect++;
181  }
182  }
183 
184  reduce(nOpen, sumOp<label>());
185  reduce(maxOpennessCell, maxOp<scalar>());
186 
187  reduce(nAspect, sumOp<label>());
188  reduce(maxAspectRatio, maxOp<scalar>());
189 
190 
191  if (nOpen > 0)
192  {
193  if (debug || report)
194  {
195  Info<< " ***Open cells found, max cell openness: "
196  << maxOpennessCell << ", number of open cells " << nOpen
197  << endl;
198  }
199 
200  return true;
201  }
202 
203  if (nAspect > 0)
204  {
205  if (debug || report)
206  {
207  Info<< " ***High aspect ratio cells found, Max aspect ratio: "
208  << maxAspectRatio
209  << ", number of cells " << nAspect
210  << endl;
211  }
212 
213  return true;
214  }
215 
216  if (debug || report)
217  {
218  Info<< " Max cell openness = " << maxOpennessCell << " OK." << nl
219  << " Max aspect ratio = " << maxAspectRatio << " OK."
220  << endl;
221  }
222 
223  return false;
224 }
225 
226 
228 (
229  const vectorField& faceAreas,
230  const bool report,
231  const bool detailedReport,
232  labelHashSet* setPtr
233 ) const
234 {
235  DebugInFunction << "Checking face area magnitudes" << endl;
236 
237  const scalarField magFaceAreas(mag(faceAreas));
238 
239  scalar minArea = GREAT;
240  scalar maxArea = -GREAT;
241 
242  forAll(magFaceAreas, facei)
243  {
244  if (magFaceAreas[facei] < VSMALL)
245  {
246  if (setPtr)
247  {
248  setPtr->insert(facei);
249  }
250  if (detailedReport)
251  {
252  if (isInternalFace(facei))
253  {
254  Pout<< "Zero or negative face area detected for "
255  << "internal face "<< facei << " between cells "
256  << faceOwner()[facei] << " and "
257  << faceNeighbour()[facei]
258  << ". Face area magnitude = " << magFaceAreas[facei]
259  << endl;
260  }
261  else
262  {
263  Pout<< "Zero or negative face area detected for "
264  << "boundary face " << facei << " next to cell "
265  << faceOwner()[facei] << ". Face area magnitude = "
266  << magFaceAreas[facei] << endl;
267  }
268  }
269  }
270 
271  minArea = min(minArea, magFaceAreas[facei]);
272  maxArea = max(maxArea, magFaceAreas[facei]);
273  }
274 
275  reduce(minArea, minOp<scalar>());
276  reduce(maxArea, maxOp<scalar>());
277 
278  if (minArea < VSMALL)
279  {
280  if (debug || report)
281  {
282  Info<< " ***Zero or negative face area detected. "
283  "Minimum area: " << minArea << endl;
284  }
285 
286  return true;
287  }
288 
289  if (debug || report)
290  {
291  Info<< " Minimum face area = " << minArea
292  << ". Maximum face area = " << maxArea
293  << ". Face area magnitudes OK." << endl;
294  }
295 
296  return false;
297 }
298 
299 
301 (
302  const scalarField& vols,
303  const bool report,
304  const bool detailedReport,
305  labelHashSet* setPtr
306 ) const
307 {
308  DebugInFunction << "Checking cell volumes" << endl;
309 
310  scalar minVolume = GREAT;
311  scalar maxVolume = -GREAT;
312 
313  label nNegVolCells = 0;
314 
315  forAll(vols, celli)
316  {
317  if (vols[celli] < VSMALL)
318  {
319  if (setPtr)
320  {
321  setPtr->insert(celli);
322  }
323  if (detailedReport)
324  {
325  Pout<< "Zero or negative cell volume detected for cell "
326  << celli << ". Volume = " << vols[celli] << endl;
327  }
328 
329  nNegVolCells++;
330  }
331 
332  minVolume = min(minVolume, vols[celli]);
333  maxVolume = max(maxVolume, vols[celli]);
334  }
335 
336  reduce(minVolume, minOp<scalar>());
337  reduce(maxVolume, maxOp<scalar>());
338  reduce(nNegVolCells, sumOp<label>());
339 
340  if (minVolume < VSMALL)
341  {
342  if (debug || report)
343  {
344  Info<< " ***Zero or negative cell volume detected. "
345  << "Minimum negative volume: " << minVolume
346  << ", Number of negative volume cells: " << nNegVolCells
347  << endl;
348  }
349 
350  return true;
351  }
352 
353  if (debug || report)
354  {
355  Info<< " Min volume = " << minVolume
356  << ". Max volume = " << maxVolume
357  << ". Total volume = " << gSum(vols)
358  << ". Cell volumes OK." << endl;
359  }
360 
361  return false;
362 }
363 
364 
366 (
367  const vectorField& fAreas,
368  const vectorField& cellCtrs,
369  const bool report,
370  labelHashSet* setPtr
371 ) const
372 {
373  DebugInFunction << "Checking mesh non-orthogonality" << endl;
374 
375  tmp<scalarField> tortho = primitiveMeshTools::faceOrthogonality
376  (
377  *this,
378  fAreas,
379  cellCtrs
380  );
381  const scalarField& ortho = tortho();
382 
383  // Severe nonorthogonality threshold
384  const scalar severeNonorthogonalityThreshold =
385  ::cos(degToRad(nonOrthThreshold_));
386 
387  scalar minDDotS = min(ortho);
388 
389  scalar sumDDotS = sum(ortho);
390 
391  label severeNonOrth = 0;
392 
393  label errorNonOrth = 0;
394 
395 
396  forAll(ortho, facei)
397  {
398  if (ortho[facei] < severeNonorthogonalityThreshold)
399  {
400  if (ortho[facei] > SMALL)
401  {
402  if (setPtr)
403  {
404  setPtr->insert(facei);
405  }
406 
407  severeNonOrth++;
408  }
409  else
410  {
411  if (setPtr)
412  {
413  setPtr->insert(facei);
414  }
415 
416  errorNonOrth++;
417  }
418  }
419  }
420 
421  reduce(minDDotS, minOp<scalar>());
422  reduce(sumDDotS, sumOp<scalar>());
423  reduce(severeNonOrth, sumOp<label>());
424  reduce(errorNonOrth, sumOp<label>());
425 
426  if (debug || report)
427  {
428  label neiSize = ortho.size();
429  reduce(neiSize, sumOp<label>());
430 
431  if (neiSize > 0)
432  {
433  if (debug || report)
434  {
435  Info<< " Mesh non-orthogonality Max: "
436  << radToDeg(::acos(minDDotS))
437  << " average: " << radToDeg(::acos(sumDDotS/neiSize))
438  << endl;
439  }
440  }
441 
442  if (severeNonOrth > 0)
443  {
444  Info<< " *Number of severely non-orthogonal faces: "
445  << severeNonOrth << "." << endl;
446  }
447  }
448 
449  if (errorNonOrth > 0)
450  {
451  if (debug || report)
452  {
453  Info<< " ***Number of non-orthogonality errors: "
454  << errorNonOrth << "." << endl;
455  }
456 
457  return true;
458  }
459 
460  if (debug || report)
461  {
462  Info<< " Non-orthogonality check OK." << endl;
463  }
464 
465  return false;
466 }
467 
468 
470 (
471  const pointField& points,
472  const vectorField& ctrs,
473  const bool report,
474  const bool detailedReport,
475  const scalar minPyrVol,
476  labelHashSet* setPtr
477 ) const
478 {
479  DebugInFunction << "Checking face orientation" << endl;
480 
481  const labelList& own = faceOwner();
482  const labelList& nei = faceNeighbour();
483  const faceList& f = faces();
484 
485 
486  scalarField ownPyrVol;
487  scalarField neiPyrVol;
488  primitiveMeshTools::facePyramidVolume
489  (
490  *this,
491  points,
492  ctrs,
493  ownPyrVol,
494  neiPyrVol
495  );
496 
497 
498  label nErrorPyrs = 0;
499 
500  forAll(ownPyrVol, facei)
501  {
502  if (ownPyrVol[facei] < minPyrVol)
503  {
504  if (setPtr)
505  {
506  setPtr->insert(facei);
507  }
508  if (detailedReport)
509  {
510  Pout<< "Negative pyramid volume: " << ownPyrVol[facei]
511  << " for face " << facei << " " << f[facei]
512  << " and owner cell: " << own[facei] << endl
513  << "Owner cell vertex labels: "
514  << cells()[own[facei]].labels(faces())
515  << endl;
516  }
517 
518  nErrorPyrs++;
519  }
520 
521  if (isInternalFace(facei))
522  {
523  if (neiPyrVol[facei] < minPyrVol)
524  {
525  if (setPtr)
526  {
527  setPtr->insert(facei);
528  }
529  if (detailedReport)
530  {
531  Pout<< "Negative pyramid volume: " << neiPyrVol[facei]
532  << " for face " << facei << " " << f[facei]
533  << " and neighbour cell: " << nei[facei] << nl
534  << "Neighbour cell vertex labels: "
535  << cells()[nei[facei]].labels(faces())
536  << endl;
537  }
538  nErrorPyrs++;
539  }
540  }
541  }
542 
543  reduce(nErrorPyrs, sumOp<label>());
544 
545  if (nErrorPyrs > 0)
546  {
547  if (debug || report)
548  {
549  Info<< " ***Error in face pyramids: "
550  << nErrorPyrs << " faces are incorrectly oriented."
551  << endl;
552  }
553 
554  return true;
555  }
556 
557  if (debug || report)
558  {
559  Info<< " Face pyramids OK." << endl;
560  }
561 
562  return false;
563 }
564 
565 
567 (
568  const pointField& points,
569  const vectorField& fCtrs,
570  const vectorField& fAreas,
571  const vectorField& cellCtrs,
572  const bool report,
573  labelHashSet* setPtr
574 ) const
575 {
576  DebugInFunction << "Checking face skewness" << endl;
577 
578  // Warn if the skew correction vector is more than skewWarning times
579  // larger than the face area vector
580 
581  tmp<scalarField> tskewness = primitiveMeshTools::faceSkewness
582  (
583  *this,
584  points,
585  fCtrs,
586  fAreas,
587  cellCtrs
588  );
589  const scalarField& skewness = tskewness();
590 
591  scalar maxSkew = max(skewness);
592  label nWarnSkew = 0;
593 
594  forAll(skewness, facei)
595  {
596  // Check if the skewness vector is greater than the PN vector.
597  // This does not cause trouble but is a good indication of a poor mesh.
598  if (skewness[facei] > skewThreshold_)
599  {
600  if (setPtr)
601  {
602  setPtr->insert(facei);
603  }
604 
605  nWarnSkew++;
606  }
607  }
608 
609  reduce(maxSkew, maxOp<scalar>());
610  reduce(nWarnSkew, sumOp<label>());
611 
612  if (nWarnSkew > 0)
613  {
614  if (debug || report)
615  {
616  Info<< " ***Max skewness = " << maxSkew
617  << ", " << nWarnSkew << " highly skew faces detected"
618  " which may impair the quality of the results"
619  << endl;
620  }
621 
622  return true;
623  }
624 
625  if (debug || report)
626  {
627  Info<< " Max skewness = " << maxSkew << " OK." << endl;
628  }
629 
630  return false;
631 }
632 
633 
635 (
636  const pointField& points,
637  const vectorField& faceAreas,
638  const bool report,
639  const scalar maxDeg,
640  labelHashSet* setPtr
641 ) const
642 {
643  DebugInFunction << "Checking face angles" << endl;
644 
645  if (maxDeg < -SMALL || maxDeg > 180+SMALL)
646  {
648  << "maxDeg should be [0..180] but is now " << maxDeg
649  << exit(FatalError);
650  }
651 
652  const scalar maxSin = Foam::sin(degToRad(maxDeg));
653 
654 
655  tmp<scalarField> tfaceAngles = primitiveMeshTools::faceConcavity
656  (
657  maxSin,
658  *this,
659  points,
660  faceAreas
661  );
662  const scalarField& faceAngles = tfaceAngles();
663 
664  scalar maxEdgeSin = max(faceAngles);
665 
666  label nConcave = 0;
667 
668  forAll(faceAngles, facei)
669  {
670  if (faceAngles[facei] > SMALL)
671  {
672  nConcave++;
673 
674  if (setPtr)
675  {
676  setPtr->insert(facei);
677  }
678  }
679  }
680 
681  reduce(nConcave, sumOp<label>());
682  reduce(maxEdgeSin, maxOp<scalar>());
683 
684  if (nConcave > 0)
685  {
686  scalar maxConcaveDegr =
687  radToDeg(Foam::asin(Foam::min(1.0, maxEdgeSin)));
688 
689  if (debug || report)
690  {
691  Info<< " *There are " << nConcave
692  << " faces with concave angles between consecutive"
693  << " edges. Max concave angle = " << maxConcaveDegr
694  << " degrees." << endl;
695  }
696 
697  return true;
698  }
699 
700  if (debug || report)
701  {
702  Info<< " All angles in faces OK." << endl;
703  }
704 
705  return false;
706 }
707 
708 
710 (
711  const pointField& points,
712  const vectorField& faceCentres,
713  const vectorField& faceAreas,
714  const bool report,
715  const scalar warnFlatness,
716  labelHashSet* setPtr
717 ) const
718 {
719  DebugInFunction << "Checking face flatness" << endl;
720 
721  if (warnFlatness < 0 || warnFlatness > 1)
722  {
724  << "warnFlatness should be [0..1] but is " << warnFlatness
725  << exit(FatalError);
726  }
727 
728  const faceList& fcs = faces();
729 
730  tmp<scalarField> tfaceFlatness = primitiveMeshTools::faceFlatness
731  (
732  *this,
733  points,
734  faceCentres,
735  faceAreas
736  );
737  const scalarField& faceFlatness = tfaceFlatness();
738 
739  scalarField magAreas(mag(faceAreas));
740 
741  scalar minFlatness = GREAT;
742  scalar sumFlatness = 0;
743  label nSummed = 0;
744  label nWarped = 0;
745 
746  forAll(faceFlatness, facei)
747  {
748  if (fcs[facei].size() > 3 && magAreas[facei] > VSMALL)
749  {
750  sumFlatness += faceFlatness[facei];
751  nSummed++;
752 
753  minFlatness = min(minFlatness, faceFlatness[facei]);
754 
755  if (faceFlatness[facei] < warnFlatness)
756  {
757  nWarped++;
758 
759  if (setPtr)
760  {
761  setPtr->insert(facei);
762  }
763  }
764  }
765  }
766 
767 
768  reduce(nWarped, sumOp<label>());
769  reduce(minFlatness, minOp<scalar>());
770 
771  reduce(nSummed, sumOp<label>());
772  reduce(sumFlatness, sumOp<scalar>());
773 
774  if (debug || report)
775  {
776  if (nSummed > 0)
777  {
778  Info<< " Face flatness (1 = flat, 0 = butterfly) : min = "
779  << minFlatness << " average = " << sumFlatness / nSummed
780  << endl;
781  }
782  }
783 
784 
785  if (nWarped > 0)
786  {
787  if (debug || report)
788  {
789  Info<< " *There are " << nWarped
790  << " faces with ratio between projected and actual area < "
791  << warnFlatness << endl;
792 
793  Info<< " Minimum ratio (minimum flatness, maximum warpage) = "
794  << minFlatness << endl;
795  }
796 
797  return true;
798  }
799 
800  if (debug || report)
801  {
802  Info<< " All face flatness OK." << endl;
803  }
804 
805  return false;
806 }
807 
808 
810 (
811  const vectorField& fAreas,
812  const pointField& fCentres,
813  const bool report,
814  labelHashSet* setPtr
815 ) const
816 {
817  DebugInFunction << "Checking for concave cells" << endl;
818 
819  const cellList& c = cells();
820  const labelList& fOwner = faceOwner();
821 
822  label nConcaveCells = 0;
823 
824  forAll(c, celli)
825  {
826  const cell& cFaces = c[celli];
827 
828  bool concave = false;
829 
830  forAll(cFaces, i)
831  {
832  if (concave)
833  {
834  break;
835  }
836 
837  label fI = cFaces[i];
838 
839  const point& fC = fCentres[fI];
840 
841  vector fN = fAreas[fI];
842 
843  fN /= max(mag(fN), VSMALL);
844 
845  // Flip normal if required so that it is always pointing out of
846  // the cell
847  if (fOwner[fI] != celli)
848  {
849  fN *= -1;
850  }
851 
852  // Is the centre of any other face of the cell on the
853  // wrong side of the plane of this face?
854 
855  forAll(cFaces, j)
856  {
857  if (j != i)
858  {
859  label fJ = cFaces[j];
860 
861  const point& pt = fCentres[fJ];
862 
863  // If the cell is concave, the point will be on the
864  // positive normal side of the plane of f, defined by
865  // its centre and normal, and the angle between (pt -
866  // fC) and fN will be less than 90 degrees, so the dot
867  // product will be positive.
868 
869  vector pC = (pt - fC);
870 
871  pC /= max(mag(pC), VSMALL);
872 
873  if ((pC & fN) > -planarCosAngle_)
874  {
875  // Concave or planar face
876 
877  concave = true;
878 
879  if (setPtr)
880  {
881  setPtr->insert(celli);
882  }
883 
884  nConcaveCells++;
885 
886  break;
887  }
888  }
889  }
890  }
891  }
892 
893  reduce(nConcaveCells, sumOp<label>());
894 
895  if (nConcaveCells > 0)
896  {
897  if (debug || report)
898  {
899  Info<< " ***Concave cells (using face planes) found,"
900  << " number of cells: " << nConcaveCells << endl;
901  }
902 
903  return true;
904  }
905 
906  if (debug || report)
907  {
908  Info<< " Concave cell check OK." << endl;
909  }
910 
911  return false;
912 }
913 
914 
916 (
917  const bool report,
918  labelHashSet* setPtr
919 ) const
920 {
921  DebugInFunction << "Checking face ordering" << endl;
922 
923  // Check whether internal faces are ordered in the upper triangular order
924  const labelList& own = faceOwner();
925  const labelList& nei = faceNeighbour();
926 
927  const cellList& c = cells();
928 
929  label internal = nInternalFaces();
930 
931  // Has error occurred?
932  bool error = false;
933  // Have multiple faces been detected?
934  label nMultipleCells = false;
935 
936  // Loop through faceCells once more and make sure that for internal cell
937  // the first label is smaller
938  for (label facei = 0; facei < internal; facei++)
939  {
940  if (own[facei] >= nei[facei])
941  {
942  error = true;
943 
944  if (setPtr)
945  {
946  setPtr->insert(facei);
947  }
948  }
949  }
950 
951  // Loop through all cells. For each cell, find the face that is internal
952  // and add it to the check list (upper triangular order).
953  // Once the list is completed, check it against the faceCell list
954 
955  forAll(c, celli)
956  {
957  const labelList& curFaces = c[celli];
958 
959  // Neighbouring cells
960  SortableList<label> nbr(curFaces.size());
961 
962  forAll(curFaces, i)
963  {
964  label facei = curFaces[i];
965 
966  if (facei >= nInternalFaces())
967  {
968  // Sort last
969  nbr[i] = labelMax;
970  }
971  else
972  {
973  label nbrCelli = nei[facei];
974 
975  if (nbrCelli == celli)
976  {
977  nbrCelli = own[facei];
978  }
979 
980  if (celli < nbrCelli)
981  {
982  // celli is master
983  nbr[i] = nbrCelli;
984  }
985  else
986  {
987  // nbrCell is master. Let it handle this face.
988  nbr[i] = labelMax;
989  }
990  }
991  }
992 
993  nbr.sort();
994 
995  // Now nbr holds the cellCells in incremental order. Check:
996  // - neighbouring cells appear only once. Since nbr is sorted this
997  // is simple check on consecutive elements
998  // - faces indexed in same order as nbr are incrementing as well.
999 
1000  label prevCell = nbr[0];
1001  label prevFace = curFaces[nbr.indices()[0]];
1002 
1003  bool hasMultipleFaces = false;
1004 
1005  for (label i = 1; i < nbr.size(); i++)
1006  {
1007  label thisCell = nbr[i];
1008  label thisFace = curFaces[nbr.indices()[i]];
1009 
1010  if (thisCell == labelMax)
1011  {
1012  break;
1013  }
1014 
1015  if (thisCell == prevCell)
1016  {
1017  hasMultipleFaces = true;
1018 
1019  if (setPtr)
1020  {
1021  setPtr->insert(prevFace);
1022  setPtr->insert(thisFace);
1023  }
1024  }
1025  else if (thisFace < prevFace)
1026  {
1027  error = true;
1028 
1029  if (setPtr)
1030  {
1031  setPtr->insert(thisFace);
1032  }
1033  }
1034 
1035  prevCell = thisCell;
1036  prevFace = thisFace;
1037  }
1038 
1039  if (hasMultipleFaces)
1040  {
1041  nMultipleCells++;
1042  }
1043  }
1044 
1045  reduce(error, orOp<bool>());
1046  reduce(nMultipleCells, sumOp<label>());
1047 
1048  if ((debug || report) && nMultipleCells > 0)
1049  {
1050  Info<< " <<Found " << nMultipleCells
1051  << " neighbouring cells with multiple inbetween faces." << endl;
1052  }
1053 
1054  if (error)
1055  {
1056  if (debug || report)
1057  {
1058  Info<< " ***Faces not in upper triangular order." << endl;
1059  }
1060 
1061  return true;
1062  }
1063 
1064  if (debug || report)
1065  {
1066  Info<< " Upper triangular ordering OK." << endl;
1067  }
1068 
1069  return false;
1070 }
1071 
1072 
1075  const bool report,
1076  labelHashSet* setPtr
1077 ) const
1078 {
1079  DebugInFunction << "Checking topological cell openness" << endl;
1080 
1081  label nOpenCells = 0;
1082 
1083  const faceList& f = faces();
1084  const cellList& c = cells();
1085 
1086  forAll(c, celli)
1087  {
1088  const labelList& curFaces = c[celli];
1089 
1090  const edgeList cellEdges = c[celli].edges(f);
1091 
1092  labelList edgeUsage(cellEdges.size(), Zero);
1093 
1094  forAll(curFaces, facei)
1095  {
1096  edgeList curFaceEdges = f[curFaces[facei]].edges();
1097 
1098  forAll(curFaceEdges, faceEdgeI)
1099  {
1100  const edge& curEdge = curFaceEdges[faceEdgeI];
1101 
1102  forAll(cellEdges, cellEdgeI)
1103  {
1104  if (cellEdges[cellEdgeI] == curEdge)
1105  {
1106  edgeUsage[cellEdgeI]++;
1107  break;
1108  }
1109  }
1110  }
1111  }
1112 
1113  edgeList singleEdges(cellEdges.size());
1114  label nSingleEdges = 0;
1115 
1116  forAll(edgeUsage, edgeI)
1117  {
1118  if (edgeUsage[edgeI] == 1)
1119  {
1120  singleEdges[nSingleEdges] = cellEdges[edgeI];
1121  nSingleEdges++;
1122  }
1123  else if (edgeUsage[edgeI] != 2)
1124  {
1125  if (setPtr)
1126  {
1127  setPtr->insert(celli);
1128  }
1129  }
1130  }
1131 
1132  if (nSingleEdges > 0)
1133  {
1134  if (setPtr)
1135  {
1136  setPtr->insert(celli);
1137  }
1138 
1139  nOpenCells++;
1140  }
1141  }
1142 
1143  reduce(nOpenCells, sumOp<label>());
1144 
1145  if (nOpenCells > 0)
1146  {
1147  if (debug || report)
1148  {
1149  Info<< " ***Open cells found, number of cells: " << nOpenCells
1150  << ". This problem may be fixable using the zipUpMesh utility."
1151  << endl;
1152  }
1153 
1154  return true;
1155  }
1156 
1157  if (debug || report)
1158  {
1159  Info<< " Topological cell zip-up check OK." << endl;
1160  }
1161 
1162  return false;
1163 }
1164 
1165 
1168  const bool report,
1169  labelHashSet* setPtr
1170 ) const
1171 {
1172  DebugInFunction << "Checking face vertices" << endl;
1173 
1174  // Check that all vertex labels are valid
1175  const faceList& f = faces();
1176 
1177  label nErrorFaces = 0;
1178 
1179  forAll(f, fI)
1180  {
1181  const face& curFace = f[fI];
1182 
1183  if (min(curFace) < 0 || max(curFace) > nPoints())
1184  {
1185  if (setPtr)
1186  {
1187  setPtr->insert(fI);
1188  }
1189 
1190  nErrorFaces++;
1191  }
1192 
1193  // Uniqueness of vertices
1194  labelHashSet facePoints(2*curFace.size());
1195 
1196  forAll(curFace, fp)
1197  {
1198  bool inserted = facePoints.insert(curFace[fp]);
1199 
1200  if (!inserted)
1201  {
1202  if (setPtr)
1203  {
1204  setPtr->insert(fI);
1205  }
1206 
1207  nErrorFaces++;
1208  }
1209  }
1210  }
1211 
1212  reduce(nErrorFaces, sumOp<label>());
1213 
1214  if (nErrorFaces > 0)
1215  {
1216  if (debug || report)
1217  {
1218  Info<< " Faces with invalid vertex labels found, "
1219  << " number of faces: " << nErrorFaces << endl;
1220  }
1221 
1222  return true;
1223  }
1224 
1225  if (debug || report)
1226  {
1227  Info<< " Face vertices OK." << endl;
1228  }
1229 
1230  return false;
1231 }
1232 
1233 
1236  const bool report,
1237  labelHashSet* setPtr
1238 ) const
1239 {
1240  DebugInFunction << "Checking points" << endl;
1241 
1242  label nFaceErrors = 0;
1243  label nCellErrors = 0;
1244 
1245  const labelListList& pf = pointFaces();
1246 
1247  forAll(pf, pointi)
1248  {
1249  if (pf[pointi].empty())
1250  {
1251  if (setPtr)
1252  {
1253  setPtr->insert(pointi);
1254  }
1255 
1256  nFaceErrors++;
1257  }
1258  }
1259 
1260 
1261  forAll(pf, pointi)
1262  {
1263  const labelList& pc = pointCells(pointi);
1264 
1265  if (pc.empty())
1266  {
1267  if (setPtr)
1268  {
1269  setPtr->insert(pointi);
1270  }
1271 
1272  nCellErrors++;
1273  }
1274  }
1275 
1276  reduce(nFaceErrors, sumOp<label>());
1277  reduce(nCellErrors, sumOp<label>());
1278 
1279  if (nFaceErrors > 0 || nCellErrors > 0)
1280  {
1281  if (debug || report)
1282  {
1283  Info<< " ***Unused points found in the mesh, "
1284  "number unused by faces: " << nFaceErrors
1285  << " number unused by cells: " << nCellErrors
1286  << endl;
1287  }
1288 
1289  return true;
1290  }
1291 
1292  if (debug || report)
1293  {
1294  Info<< " Point usage OK." << endl;
1295  }
1296 
1297  return false;
1298 }
1299 
1300 
1303  const label facei,
1304  const Map<label>& nCommonPoints,
1305  label& nBaffleFaces,
1306  labelHashSet* setPtr
1307 ) const
1308 {
1309  bool error = false;
1310 
1311  forAllConstIters(nCommonPoints, iter)
1312  {
1313  label nbFacei = iter.key();
1314  label nCommon = iter();
1315 
1316  const face& curFace = faces()[facei];
1317  const face& nbFace = faces()[nbFacei];
1318 
1319  if (nCommon == nbFace.size() || nCommon == curFace.size())
1320  {
1321  if (nbFace.size() != curFace.size())
1322  {
1323  error = true;
1324  }
1325  else
1326  {
1327  nBaffleFaces++;
1328  }
1329 
1330  if (setPtr)
1331  {
1332  setPtr->insert(facei);
1333  setPtr->insert(nbFacei);
1334  }
1335  }
1336  }
1337 
1338  return error;
1339 }
1340 
1341 
1344  const label facei,
1345  const Map<label>& nCommonPoints,
1346  labelHashSet* setPtr
1347 ) const
1348 {
1349  bool error = false;
1350 
1351  forAllConstIters(nCommonPoints, iter)
1352  {
1353  label nbFacei = iter.key();
1354  label nCommon = iter();
1355 
1356  const face& curFace = faces()[facei];
1357  const face& nbFace = faces()[nbFacei];
1358 
1359  if
1360  (
1361  nCommon >= 2
1362  && nCommon != nbFace.size()
1363  && nCommon != curFace.size()
1364  )
1365  {
1366  forAll(curFace, fp)
1367  {
1368  // Get the index in the neighbouring face shared with curFace
1369  label nb = nbFace.find(curFace[fp]);
1370 
1371  if (nb != -1)
1372  {
1373 
1374  // Check the whole face from nb onwards for shared vertices
1375  // with neighbouring face. Rule is that any shared vertices
1376  // should be consecutive on both faces i.e. if they are
1377  // vertices fp,fp+1,fp+2 on one face they should be
1378  // vertices nb, nb+1, nb+2 (or nb+2, nb+1, nb) on the
1379  // other face.
1380 
1381 
1382  // Vertices before and after on curFace
1383  label fpPlus1 = curFace.fcIndex(fp);
1384  label fpMin1 = curFace.rcIndex(fp);
1385 
1386  // Vertices before and after on nbFace
1387  label nbPlus1 = nbFace.fcIndex(nb);
1388  label nbMin1 = nbFace.rcIndex(nb);
1389 
1390  // Find order of walking by comparing next points on both
1391  // faces.
1392  label curInc = labelMax;
1393  label nbInc = labelMax;
1394 
1395  if (nbFace[nbPlus1] == curFace[fpPlus1])
1396  {
1397  curInc = 1;
1398  nbInc = 1;
1399  }
1400  else if (nbFace[nbPlus1] == curFace[fpMin1])
1401  {
1402  curInc = -1;
1403  nbInc = 1;
1404  }
1405  else if (nbFace[nbMin1] == curFace[fpMin1])
1406  {
1407  curInc = -1;
1408  nbInc = -1;
1409  }
1410  else
1411  {
1412  curInc = 1;
1413  nbInc = -1;
1414  }
1415 
1416 
1417  // Pass1: loop until start of common vertices found.
1418  label curNb = nb;
1419  label curFp = fp;
1420 
1421  do
1422  {
1423  curFp += curInc;
1424 
1425  if (curFp >= curFace.size())
1426  {
1427  curFp = 0;
1428  }
1429  else if (curFp < 0)
1430  {
1431  curFp = curFace.size()-1;
1432  }
1433 
1434  curNb += nbInc;
1435 
1436  if (curNb >= nbFace.size())
1437  {
1438  curNb = 0;
1439  }
1440  else if (curNb < 0)
1441  {
1442  curNb = nbFace.size()-1;
1443  }
1444  } while (curFace[curFp] == nbFace[curNb]);
1445 
1446 
1447  // Pass2: check equality walking from curFp, curNb
1448  // in opposite order.
1449 
1450  curInc = -curInc;
1451  nbInc = -nbInc;
1452 
1453  for (label commonI = 0; commonI < nCommon; commonI++)
1454  {
1455  curFp += curInc;
1456 
1457  if (curFp >= curFace.size())
1458  {
1459  curFp = 0;
1460  }
1461  else if (curFp < 0)
1462  {
1463  curFp = curFace.size()-1;
1464  }
1465 
1466  curNb += nbInc;
1467 
1468  if (curNb >= nbFace.size())
1469  {
1470  curNb = 0;
1471  }
1472  else if (curNb < 0)
1473  {
1474  curNb = nbFace.size()-1;
1475  }
1476 
1477  if (curFace[curFp] != nbFace[curNb])
1478  {
1479  if (setPtr)
1480  {
1481  setPtr->insert(facei);
1482  setPtr->insert(nbFacei);
1483  }
1484 
1485  error = true;
1486 
1487  break;
1488  }
1489  }
1490 
1491 
1492  // Done the curFace - nbFace combination.
1493  break;
1494  }
1495  }
1496  }
1497  }
1498 
1499  return error;
1500 }
1501 
1502 
1505  const bool report,
1506  labelHashSet* setPtr
1507 ) const
1508 {
1509  DebugInFunction << "Checking face-face connectivity" << endl;
1510 
1511  const labelListList& pf = pointFaces();
1512 
1513  label nBaffleFaces = 0;
1514  label nErrorDuplicate = 0;
1515  label nErrorOrder = 0;
1516  Map<label> nCommonPoints(128);
1517 
1518  for (label facei = 0; facei < nFaces(); facei++)
1519  {
1520  const face& curFace = faces()[facei];
1521 
1522  // Calculate number of common points between current facei and
1523  // neighbouring face. Store on map.
1524  nCommonPoints.clear();
1525 
1526  forAll(curFace, fp)
1527  {
1528  label pointi = curFace[fp];
1529 
1530  const labelList& nbs = pf[pointi];
1531 
1532  forAll(nbs, nbI)
1533  {
1534  label nbFacei = nbs[nbI];
1535 
1536  if (facei < nbFacei)
1537  {
1538  // Only check once for each combination of two faces.
1539  ++(nCommonPoints(nbFacei, 0));
1540  }
1541  }
1542  }
1543 
1544  // Perform various checks on common points
1545 
1546  // Check all vertices shared (duplicate point)
1547  if (checkDuplicateFaces(facei, nCommonPoints, nBaffleFaces, setPtr))
1548  {
1549  nErrorDuplicate++;
1550  }
1551 
1552  // Check common vertices are consecutive on both faces
1553  if (checkCommonOrder(facei, nCommonPoints, setPtr))
1554  {
1555  nErrorOrder++;
1556  }
1557  }
1558 
1559  reduce(nBaffleFaces, sumOp<label>());
1560  reduce(nErrorDuplicate, sumOp<label>());
1561  reduce(nErrorOrder, sumOp<label>());
1562 
1563  if (nBaffleFaces)
1564  {
1565  Info<< " Number of identical duplicate faces (baffle faces): "
1566  << nBaffleFaces << endl;
1567  }
1568 
1569  if (nErrorDuplicate > 0 || nErrorOrder > 0)
1570  {
1571  // These are actually warnings, not errors.
1572  if (nErrorDuplicate > 0)
1573  {
1574  Info<< " <<Number of duplicate (not baffle) faces found: "
1575  << nErrorDuplicate
1576  << ". This might indicate a problem." << endl;
1577  }
1578 
1579  if (nErrorOrder > 0)
1580  {
1581  Info<< " <<Number of faces with non-consecutive shared points: "
1582  << nErrorOrder << ". This might indicate a problem." << endl;
1583  }
1584 
1585  return false; //return true;
1586  }
1587 
1588  if (debug || report)
1589  {
1590  Info<< " Face-face connectivity OK." << endl;
1591  }
1592 
1593  return false;
1594 }
1595 
1596 
1597 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1598 
1599 bool Foam::primitiveMesh::checkClosedBoundary(const bool report) const
1600 {
1601  return checkClosedBoundary(faceAreas(), report, bitSet());
1602 }
1603 
1604 
1607  const bool report,
1608  labelHashSet* setPtr,
1609  labelHashSet* aspectSetPtr,
1610  const Vector<label>& solutionD
1611 ) const
1612 {
1613  return checkClosedCells
1614  (
1615  faceAreas(),
1616  cellVolumes(),
1617  report,
1618  setPtr,
1619  aspectSetPtr,
1620  solutionD
1621  );
1622 }
1623 
1624 
1627  const bool report,
1628  labelHashSet* setPtr
1629 ) const
1630 {
1631  return checkFaceAreas
1632  (
1633  faceAreas(),
1634  report,
1635  false, // detailedReport,
1636  setPtr
1637  );
1638 }
1639 
1640 
1643  const bool report,
1644  labelHashSet* setPtr
1645 ) const
1646 {
1647  return checkCellVolumes
1648  (
1649  cellVolumes(),
1650  report,
1651  false, // detailedReport,
1652  setPtr
1653  );
1654 }
1655 
1656 
1659  const bool report,
1660  labelHashSet* setPtr
1661 ) const
1662 {
1663  return checkFaceOrthogonality
1664  (
1665  faceAreas(),
1666  cellCentres(),
1667  report,
1668  setPtr
1669  );
1670 }
1671 
1672 
1675  const bool report,
1676  const scalar minPyrVol,
1677  labelHashSet* setPtr
1678 ) const
1679 {
1680  return checkFacePyramids
1681  (
1682  points(),
1683  cellCentres(),
1684  report,
1685  false, // detailedReport,
1686  minPyrVol,
1687  setPtr
1688  );
1689 }
1690 
1691 
1694  const bool report,
1695  labelHashSet* setPtr
1696 ) const
1697 {
1698  return checkFaceSkewness
1699  (
1700  points(),
1701  faceCentres(),
1702  faceAreas(),
1703  cellCentres(),
1704  report,
1705  setPtr
1706  );
1707 }
1708 
1709 
1712  const bool report,
1713  const scalar maxDeg,
1714  labelHashSet* setPtr
1715 ) const
1716 {
1717  return checkFaceAngles
1718  (
1719  points(),
1720  faceAreas(),
1721  report,
1722  maxDeg,
1723  setPtr
1724  );
1725 }
1726 
1727 
1730  const bool report,
1731  const scalar warnFlatness,
1732  labelHashSet* setPtr
1733 ) const
1734 {
1735  return checkFaceFlatness
1736  (
1737  points(),
1738  faceCentres(),
1739  faceAreas(),
1740  report,
1741  warnFlatness,
1742  setPtr
1743  );
1744 }
1745 
1746 
1749  const bool report,
1750  labelHashSet* setPtr
1751 ) const
1752 {
1753  return checkConcaveCells
1754  (
1755  faceAreas(),
1756  faceCentres(),
1757  report,
1758  setPtr
1759  );
1760 }
1761 
1762 
1763 bool Foam::primitiveMesh::checkTopology(const bool report) const
1764 {
1765  label nFailedChecks = 0;
1766 
1767  if (checkPoints(report)) ++nFailedChecks;
1768  if (checkUpperTriangular(report)) ++nFailedChecks;
1769  if (checkCellsZipUp(report)) ++nFailedChecks;
1770  if (checkFaceVertices(report)) ++nFailedChecks;
1771  if (checkFaceFaces(report)) ++nFailedChecks;
1772 
1773  if (nFailedChecks)
1774  {
1775  if (debug || report)
1776  {
1777  Info<< " Failed " << nFailedChecks
1778  << " mesh topology checks." << endl;
1779  }
1780 
1781  return true;
1782  }
1783 
1784  if (debug || report)
1785  {
1786  Info<< " Mesh topology OK." << endl;
1787  }
1788 
1789  return false;
1790 }
1791 
1792 
1793 bool Foam::primitiveMesh::checkGeometry(const bool report) const
1794 {
1795  label nFailedChecks = 0;
1796 
1797  if (checkClosedBoundary(report)) ++nFailedChecks;
1798  if (checkClosedCells(report)) ++nFailedChecks;
1799  if (checkFaceAreas(report)) ++nFailedChecks;
1800  if (checkCellVolumes(report)) ++nFailedChecks;
1801  if (checkFaceOrthogonality(report)) ++nFailedChecks;
1802  if (checkFacePyramids(report)) ++nFailedChecks;
1803  if (checkFaceSkewness(report)) ++nFailedChecks;
1804 
1805  if (nFailedChecks)
1806  {
1807  if (debug || report)
1808  {
1809  Info<< " Failed " << nFailedChecks
1810  << " mesh geometry checks." << endl;
1811  }
1812 
1813  return true;
1814  }
1815 
1816  if (debug || report)
1817  {
1818  Info<< " Mesh geometry OK." << endl;
1819  }
1820 
1821  return false;
1822 }
1823 
1824 
1825 bool Foam::primitiveMesh::checkMesh(const bool report) const
1826 {
1827  DebugInFunction << "Checking primitiveMesh" << endl;
1828 
1829  label nFailedChecks = checkTopology(report) + checkGeometry(report);
1830 
1831  if (nFailedChecks)
1832  {
1833  if (debug || report)
1834  {
1835  Info<< " Failed " << nFailedChecks
1836  << " mesh checks." << endl;
1837  }
1838 
1839  return true;
1840  }
1841 
1842  if (debug || report)
1843  {
1844  Info<< "Mesh OK." << endl;
1845  }
1846 
1847  return false;
1848 }
1849 
1850 
1851 Foam::scalar Foam::primitiveMesh::setClosedThreshold(const scalar val)
1852 {
1853  scalar prev = closedThreshold_;
1854  closedThreshold_ = val;
1855 
1856  return prev;
1857 }
1858 
1859 
1860 Foam::scalar Foam::primitiveMesh::setAspectThreshold(const scalar val)
1861 {
1862  scalar prev = aspectThreshold_;
1863  aspectThreshold_ = val;
1864 
1865  return prev;
1866 }
1867 
1868 
1869 Foam::scalar Foam::primitiveMesh::setNonOrthThreshold(const scalar val)
1870 {
1871  scalar prev = nonOrthThreshold_;
1872  nonOrthThreshold_ = val;
1873 
1874  return prev;
1875 }
1876 
1877 
1878 Foam::scalar Foam::primitiveMesh::setSkewThreshold(const scalar val)
1879 {
1880  scalar prev = skewThreshold_;
1881  skewThreshold_ = val;
1882 
1883  return prev;
1884 }
1885 
1886 
1887 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::primitiveMesh::checkFaceFlatness
bool checkFaceFlatness(const pointField &points, const vectorField &faceCentres, const vectorField &faceAreas, const bool report, const scalar warnFlatness, labelHashSet *setPtr) const
Check face warpage.
Definition: primitiveMeshCheck.C:710
Foam::primitiveMesh::checkFaceSkewness
bool checkFaceSkewness(const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs, const bool report, labelHashSet *setPtr) const
Check face skewness.
Definition: primitiveMeshCheck.C:567
Foam::maxOp
Definition: ops.H:223
Foam::primitiveMesh::checkFacePyramids
bool checkFacePyramids(const pointField &points, const vectorField &ctrs, const bool report, const bool detailedReport, const scalar minPyrVol, labelHashSet *setPtr) const
Check face pyramid volume.
Definition: primitiveMeshCheck.C:470
Foam::primitiveMesh::checkFaceAreas
bool checkFaceAreas(const vectorField &faceAreas, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative face areas.
Definition: primitiveMeshCheck.C:228
Foam::primitiveMesh::checkFaceFaces
virtual bool checkFaceFaces(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face-face connectivity.
Definition: primitiveMeshCheck.C:1504
Foam::labelMax
constexpr label labelMax
Definition: label.H:61
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::primitiveMesh::checkMesh
virtual bool checkMesh(const bool report=false) const
Check mesh for correctness. Returns false for no error.
Definition: primitiveMeshCheck.C:1825
primitiveMeshTools.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::minOp
Definition: ops.H:224
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::primitiveMesh::setSkewThreshold
static scalar setSkewThreshold(const scalar)
Set the skewness warning threshold as percentage.
Definition: primitiveMeshCheck.C:1878
Foam::primitiveMesh::checkDuplicateFaces
bool checkDuplicateFaces(const label, const Map< label > &, label &nBaffleFaces, labelHashSet *) const
Check if all points on face are shared with another face.
Definition: primitiveMeshCheck.C:1302
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::primitiveMesh::checkFaceAngles
bool checkFaceAngles(const pointField &points, const vectorField &faceAreas, const bool report, const scalar maxDeg, labelHashSet *setPtr) const
Check face angles.
Definition: primitiveMeshCheck.C:635
Foam::Map< label >
unitConversion.H
Unit conversion functions.
Foam::primitiveMesh::checkCellVolumes
bool checkCellVolumes(const scalarField &vols, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative cell volumes.
Definition: primitiveMeshCheck.C:301
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::checkTopology
label checkTopology(const polyMesh &mesh, const bool allTopology, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, const autoPtr< writer< scalar >> &setWriter)
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::HashSet< label, Hash< label > >
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::primitiveMesh::checkTopology
virtual bool checkTopology(const bool report=false) const
Check mesh topology for correctness.
Definition: primitiveMeshCheck.C:1763
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::primitiveMesh::checkConcaveCells
bool checkConcaveCells(const vectorField &fAreas, const pointField &fCentres, const bool report, labelHashSet *setPtr) const
Check for concave cells by the planes of faces.
Definition: primitiveMeshCheck.C:810
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
Foam::primitiveMesh::setClosedThreshold
static scalar setClosedThreshold(const scalar)
Set the closedness ratio warning threshold.
Definition: primitiveMeshCheck.C:1851
Foam::cmptMax
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:253
SortableList.H
Foam::primitiveMesh::checkGeometry
virtual bool checkGeometry(const bool report=false) const
Check mesh geometry (& implicitly topology) for correctness.
Definition: primitiveMeshCheck.C:1793
Foam::primitiveMesh::setNonOrthThreshold
static scalar setNonOrthThreshold(const scalar)
Set the non-orthogonality warning threshold in degrees.
Definition: primitiveMeshCheck.C:1869
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::radToDeg
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
Definition: unitConversion.H:54
Foam::checkGeometry
label checkGeometry(const polyMesh &mesh, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, const autoPtr< writer< scalar >> &setWriter)
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::FatalError
error FatalError
reduce
reduce(hasMovingMesh, orOp< bool >())
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:60
Foam::primitiveMesh::aspectThreshold_
static scalar aspectThreshold_
Aspect ratio warning threshold.
Definition: primitiveMesh.H:241
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::primitiveMesh::checkPoints
virtual bool checkPoints(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for unused points.
Definition: primitiveMeshCheck.C:1235
Foam::primitiveMesh::closedThreshold_
static scalar closedThreshold_
Static data to control mesh checking.
Definition: primitiveMesh.H:238
Foam::primitiveMesh::nonOrthThreshold_
static scalar nonOrthThreshold_
Non-orthogonality warning threshold in deg.
Definition: primitiveMesh.H:244
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
edgeHashes.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::primitiveMesh::planarCosAngle_
static scalar planarCosAngle_
Threshold where faces are considered coplanar.
Definition: primitiveMesh.H:250
Foam::primitiveMesh::setAspectThreshold
static scalar setAspectThreshold(const scalar)
Set the aspect ratio warning threshold.
Definition: primitiveMeshCheck.C:1860
f
labelList f(nPoints)
Foam::primitiveMesh::skewThreshold_
static scalar skewThreshold_
Skewness warning threshold.
Definition: primitiveMesh.H:247
Foam::primitiveMesh::checkClosedCells
bool checkClosedCells(const vectorField &faceAreas, const scalarField &cellVolumes, const bool report, labelHashSet *setPtr, labelHashSet *aspectSetPtr, const Vector< label > &meshD) const
Check cells for closedness.
Definition: primitiveMeshCheck.C:100
Foam::Vector
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:62
Foam::primitiveMesh::checkClosedBoundary
bool checkClosedBoundary(const vectorField &areas, const bool report, const bitSet &internalOrCoupledFaces) const
Check boundary for closedness.
Definition: primitiveMeshCheck.C:49
Foam::List< cell >
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::PackedList::size
label size() const noexcept
Number of entries.
Definition: PackedListI.H:377
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::primitiveMesh::checkFaceOrthogonality
bool checkFaceOrthogonality(const vectorField &fAreas, const vectorField &cellCtrs, const bool report, labelHashSet *setPtr) const
Check for non-orthogonality.
Definition: primitiveMeshCheck.C:366
Foam::pointCells
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition: pointCells.H:56
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
ListOps.H
Various functions to operate on Lists.
Foam::orOp
Definition: ops.H:234
Foam::primitiveMesh::checkFaceVertices
virtual bool checkFaceVertices(const bool report=false, labelHashSet *setPtr=nullptr) const
Check uniqueness of face vertices.
Definition: primitiveMeshCheck.C:1167
pyramidPointFaceRef.H
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::error
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:73
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:267
Foam::primitiveMesh::checkCellsZipUp
virtual bool checkCellsZipUp(const bool report=false, labelHashSet *setPtr=nullptr) const
Check cell zip-up.
Definition: primitiveMeshCheck.C:1074
Foam::primitiveMesh::checkCommonOrder
bool checkCommonOrder(const label, const Map< label > &, labelHashSet *) const
Check that shared points are in consecutive order.
Definition: primitiveMeshCheck.C:1343
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:89
Foam::primitiveMesh::checkUpperTriangular
virtual bool checkUpperTriangular(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face ordering.
Definition: primitiveMeshCheck.C:916