Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "pyramid.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
39Foam::scalar Foam::primitiveMesh::closedThreshold_ = 1.0e-6;
43Foam::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 solveVector 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 += areas[facei];
68 sumMagClosedBoundary += mag(areas[facei]);
69 }
70 }
71
72 reduce(sumClosed, sumOp<solveVector>());
73 reduce(sumMagClosedBoundary, sumOp<solveScalar>());
74
75 solveVector 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;
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."
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 const label neiSize = returnReduce(ortho.size(), sumOp<label>());
429
430 if (neiSize > 0)
431 {
432 if (debug || report)
433 {
434 Info<< " Mesh non-orthogonality Max: "
435 << radToDeg(::acos(minDDotS))
436 << " average: " << radToDeg(::acos(sumDDotS/neiSize))
437 << endl;
438 }
439 }
440
441 if (severeNonOrth > 0)
442 {
443 Info<< " *Number of severely non-orthogonal faces: "
444 << severeNonOrth << "." << endl;
445 }
446 }
447
448 if (errorNonOrth > 0)
449 {
450 if (debug || report)
451 {
452 Info<< " ***Number of non-orthogonality errors: "
453 << errorNonOrth << "." << endl;
454 }
455
456 return true;
457 }
458
459 if (debug || report)
460 {
461 Info<< " Non-orthogonality check OK." << endl;
462 }
463
464 return false;
465}
466
467
469(
470 const pointField& points,
471 const vectorField& ctrs,
472 const bool report,
473 const bool detailedReport,
474 const scalar minPyrVol,
475 labelHashSet* setPtr
476) const
477{
478 DebugInFunction << "Checking face orientation" << endl;
479
480 const labelList& own = faceOwner();
481 const labelList& nei = faceNeighbour();
482 const faceList& f = faces();
483
484
485 scalarField ownPyrVol;
486 scalarField neiPyrVol;
488 (
489 *this,
490 points,
491 ctrs,
492 ownPyrVol,
493 neiPyrVol
494 );
495
496
497 label nErrorPyrs = 0;
498
499 forAll(ownPyrVol, facei)
500 {
501 if (ownPyrVol[facei] < minPyrVol)
502 {
503 if (setPtr)
504 {
505 setPtr->insert(facei);
506 }
507 if (detailedReport)
508 {
509 Pout<< "Negative pyramid volume: " << ownPyrVol[facei]
510 << " for face " << facei << " " << f[facei]
511 << " and owner cell: " << own[facei] << endl
512 << "Owner cell vertex labels: "
513 << cells()[own[facei]].labels(faces())
514 << endl;
515 }
516
517 nErrorPyrs++;
518 }
519
520 if (isInternalFace(facei))
521 {
522 if (neiPyrVol[facei] < minPyrVol)
523 {
524 if (setPtr)
525 {
526 setPtr->insert(facei);
527 }
528 if (detailedReport)
529 {
530 Pout<< "Negative pyramid volume: " << neiPyrVol[facei]
531 << " for face " << facei << " " << f[facei]
532 << " and neighbour cell: " << nei[facei] << nl
533 << "Neighbour cell vertex labels: "
534 << cells()[nei[facei]].labels(faces())
535 << endl;
536 }
537 nErrorPyrs++;
538 }
539 }
540 }
541
542 reduce(nErrorPyrs, sumOp<label>());
543
544 if (nErrorPyrs > 0)
545 {
546 if (debug || report)
547 {
548 Info<< " ***Error in face pyramids: "
549 << nErrorPyrs << " faces are incorrectly oriented."
550 << endl;
551 }
552
553 return true;
554 }
555
556 if (debug || report)
557 {
558 Info<< " Face pyramids OK." << endl;
559 }
560
561 return false;
562}
563
564
566(
567 const pointField& points,
568 const vectorField& fCtrs,
569 const vectorField& fAreas,
570 const vectorField& cellCtrs,
571 const bool report,
572 labelHashSet* setPtr
573) const
574{
575 DebugInFunction << "Checking face skewness" << endl;
576
577 // Warn if the skew correction vector is more than skewWarning times
578 // larger than the face area vector
579
580 tmp<scalarField> tskewness = primitiveMeshTools::faceSkewness
581 (
582 *this,
583 points,
584 fCtrs,
585 fAreas,
586 cellCtrs
587 );
588 const scalarField& skewness = tskewness();
589
590 scalar maxSkew = max(skewness);
591 label nWarnSkew = 0;
592
593 forAll(skewness, facei)
594 {
595 // Check if the skewness vector is greater than the PN vector.
596 // This does not cause trouble but is a good indication of a poor mesh.
597 if (skewness[facei] > skewThreshold_)
598 {
599 if (setPtr)
600 {
601 setPtr->insert(facei);
602 }
603
604 nWarnSkew++;
605 }
606 }
607
608 reduce(maxSkew, maxOp<scalar>());
609 reduce(nWarnSkew, sumOp<label>());
610
611 if (nWarnSkew > 0)
612 {
613 if (debug || report)
614 {
615 Info<< " ***Max skewness = " << maxSkew
616 << ", " << nWarnSkew << " highly skew faces detected"
617 " which may impair the quality of the results"
618 << endl;
619 }
620
621 return true;
622 }
623
624 if (debug || report)
625 {
626 Info<< " Max skewness = " << maxSkew << " OK." << endl;
627 }
628
629 return false;
630}
631
632
634(
635 const pointField& points,
636 const vectorField& faceAreas,
637 const bool report,
638 const scalar maxDeg,
639 labelHashSet* setPtr
640) const
641{
642 DebugInFunction << "Checking face angles" << endl;
643
644 if (maxDeg < -SMALL || maxDeg > 180+SMALL)
645 {
647 << "maxDeg should be [0..180] but is now " << maxDeg
648 << exit(FatalError);
649 }
650
651 const scalar maxSin = Foam::sin(degToRad(maxDeg));
652
653
655 (
656 maxSin,
657 *this,
658 points,
659 faceAreas
660 );
661 const scalarField& faceAngles = tfaceAngles();
662
663 scalar maxEdgeSin = max(faceAngles);
664
665 label nConcave = 0;
666
667 forAll(faceAngles, facei)
668 {
669 if (faceAngles[facei] > SMALL)
670 {
671 nConcave++;
672
673 if (setPtr)
674 {
675 setPtr->insert(facei);
676 }
677 }
678 }
679
680 reduce(nConcave, sumOp<label>());
681 reduce(maxEdgeSin, maxOp<scalar>());
682
683 if (nConcave > 0)
684 {
685 scalar maxConcaveDegr =
686 radToDeg(Foam::asin(Foam::min(1.0, maxEdgeSin)));
687
688 if (debug || report)
689 {
690 Info<< " *There are " << nConcave
691 << " faces with concave angles between consecutive"
692 << " edges. Max concave angle = " << maxConcaveDegr
693 << " degrees." << endl;
694 }
695
696 return true;
697 }
698
699 if (debug || report)
700 {
701 Info<< " All angles in faces OK." << endl;
702 }
703
704 return false;
705}
706
707
709(
710 const pointField& points,
711 const vectorField& faceCentres,
712 const vectorField& faceAreas,
713 const bool report,
714 const scalar warnFlatness,
715 labelHashSet* setPtr
716) const
717{
718 DebugInFunction << "Checking face flatness" << endl;
719
720 if (warnFlatness < 0 || warnFlatness > 1)
721 {
723 << "warnFlatness should be [0..1] but is " << warnFlatness
724 << exit(FatalError);
725 }
726
727 const faceList& fcs = faces();
728
730 (
731 *this,
732 points,
733 faceCentres,
734 faceAreas
735 );
736 const scalarField& faceFlatness = tfaceFlatness();
737
738 scalarField magAreas(mag(faceAreas));
739
740 scalar minFlatness = GREAT;
741 scalar sumFlatness = 0;
742 label nSummed = 0;
743 label nWarped = 0;
744
745 forAll(faceFlatness, facei)
746 {
747 if (fcs[facei].size() > 3 && magAreas[facei] > VSMALL)
748 {
749 sumFlatness += faceFlatness[facei];
750 nSummed++;
751
752 minFlatness = min(minFlatness, faceFlatness[facei]);
753
754 if (faceFlatness[facei] < warnFlatness)
755 {
756 nWarped++;
757
758 if (setPtr)
759 {
760 setPtr->insert(facei);
761 }
762 }
763 }
764 }
765
766
767 reduce(nWarped, sumOp<label>());
768 reduce(minFlatness, minOp<scalar>());
769
770 reduce(nSummed, sumOp<label>());
771 reduce(sumFlatness, sumOp<scalar>());
772
773 if (debug || report)
774 {
775 if (nSummed > 0)
776 {
777 Info<< " Face flatness (1 = flat, 0 = butterfly) : min = "
778 << minFlatness << " average = " << sumFlatness / nSummed
779 << endl;
780 }
781 }
782
783
784 if (nWarped > 0)
785 {
786 if (debug || report)
787 {
788 Info<< " *There are " << nWarped
789 << " faces with ratio between projected and actual area < "
790 << warnFlatness << endl;
791
792 Info<< " Minimum ratio (minimum flatness, maximum warpage) = "
793 << minFlatness << endl;
794 }
795
796 return true;
797 }
798
799 if (debug || report)
800 {
801 Info<< " All face flatness OK." << endl;
802 }
803
804 return false;
805}
806
807
809(
810 const vectorField& fAreas,
811 const pointField& fCentres,
812 const bool report,
813 labelHashSet* setPtr
814) const
815{
816 DebugInFunction << "Checking for concave cells" << endl;
817
818 const cellList& c = cells();
819 const labelList& fOwner = faceOwner();
820
821 label nConcaveCells = 0;
822
823 forAll(c, celli)
824 {
825 const cell& cFaces = c[celli];
826
827 bool concave = false;
828
829 forAll(cFaces, i)
830 {
831 if (concave)
832 {
833 break;
834 }
835
836 label fI = cFaces[i];
837
838 const point& fC = fCentres[fI];
839
840 vector fN = fAreas[fI];
841
842 fN /= max(mag(fN), VSMALL);
843
844 // Flip normal if required so that it is always pointing out of
845 // the cell
846 if (fOwner[fI] != celli)
847 {
848 fN *= -1;
849 }
850
851 // Is the centre of any other face of the cell on the
852 // wrong side of the plane of this face?
853
854 forAll(cFaces, j)
855 {
856 if (j != i)
857 {
858 label fJ = cFaces[j];
859
860 const point& pt = fCentres[fJ];
861
862 // If the cell is concave, the point will be on the
863 // positive normal side of the plane of f, defined by
864 // its centre and normal, and the angle between (pt -
865 // fC) and fN will be less than 90 degrees, so the dot
866 // product will be positive.
867
868 vector pC = (pt - fC);
869
870 pC /= max(mag(pC), VSMALL);
871
872 if ((pC & fN) > -planarCosAngle_)
873 {
874 // Concave or planar face
875
876 concave = true;
877
878 if (setPtr)
879 {
880 setPtr->insert(celli);
881 }
882
883 nConcaveCells++;
884
885 break;
886 }
887 }
888 }
889 }
890 }
891
892 reduce(nConcaveCells, sumOp<label>());
893
894 if (nConcaveCells > 0)
895 {
896 if (debug || report)
897 {
898 Info<< " ***Concave cells (using face planes) found,"
899 << " number of cells: " << nConcaveCells << endl;
900 }
901
902 return true;
903 }
904
905 if (debug || report)
906 {
907 Info<< " Concave cell check OK." << endl;
908 }
909
910 return false;
911}
912
913
915(
916 const bool report,
917 labelHashSet* setPtr
918) const
919{
920 DebugInFunction << "Checking face ordering" << endl;
921
922 // Check whether internal faces are ordered in the upper triangular order
923 const labelList& own = faceOwner();
924 const labelList& nei = faceNeighbour();
925
926 const cellList& c = cells();
927
928 label internal = nInternalFaces();
929
930 // Has error occurred?
931 bool hasError = false;
932 // Multiple faces detected?
933 label nMultipleCells = 0;
934
935 // Loop through faceCells once more and make sure that for internal cell
936 // the first label is smaller
937 for (label facei = 0; facei < internal; facei++)
938 {
939 if (own[facei] >= nei[facei])
940 {
941 hasError = true;
942
943 if (setPtr)
944 {
945 setPtr->insert(facei);
946 }
947 }
948 }
949
950 // Loop through all cells. For each cell, find the face that is internal
951 // and add it to the check list (upper triangular order).
952 // Once the list is completed, check it against the faceCell list
953
954 forAll(c, celli)
955 {
956 const labelList& curFaces = c[celli];
957
958 // Neighbouring cells
959 SortableList<label> nbr(curFaces.size());
960
961 forAll(curFaces, i)
962 {
963 label facei = curFaces[i];
964
965 if (facei >= nInternalFaces())
966 {
967 // Sort last
968 nbr[i] = labelMax;
969 }
970 else
971 {
972 label nbrCelli = nei[facei];
973
974 if (nbrCelli == celli)
975 {
976 nbrCelli = own[facei];
977 }
978
979 if (celli < nbrCelli)
980 {
981 // celli is master
982 nbr[i] = nbrCelli;
983 }
984 else
985 {
986 // nbrCell is master. Let it handle this face.
987 nbr[i] = labelMax;
988 }
989 }
990 }
991
992 nbr.sort();
993
994 // Now nbr holds the cellCells in incremental order. Check:
995 // - neighbouring cells appear only once. Since nbr is sorted this
996 // is simple check on consecutive elements
997 // - faces indexed in same order as nbr are incrementing as well.
998
999 label prevCell = nbr[0];
1000 label prevFace = curFaces[nbr.indices()[0]];
1001
1002 bool hasMultipleFaces = false;
1003
1004 for (label i = 1; i < nbr.size(); i++)
1005 {
1006 label thisCell = nbr[i];
1007 label thisFace = curFaces[nbr.indices()[i]];
1008
1009 if (thisCell == labelMax)
1010 {
1011 break;
1012 }
1013
1014 if (thisCell == prevCell)
1015 {
1016 hasMultipleFaces = true;
1017
1018 if (setPtr)
1019 {
1020 setPtr->insert(prevFace);
1021 setPtr->insert(thisFace);
1022 }
1023 }
1024 else if (thisFace < prevFace)
1025 {
1026 hasError = true;
1027
1028 if (setPtr)
1029 {
1030 setPtr->insert(thisFace);
1031 }
1032 }
1033
1034 prevCell = thisCell;
1035 prevFace = thisFace;
1036 }
1037
1038 if (hasMultipleFaces)
1039 {
1040 nMultipleCells++;
1041 }
1042 }
1043
1044 Pstream::reduceOr(hasError);
1045 reduce(nMultipleCells, sumOp<label>());
1046
1047 if ((debug || report) && nMultipleCells > 0)
1048 {
1049 Info<< " <<Found " << nMultipleCells
1050 << " neighbouring cells with multiple inbetween faces." << endl;
1051 }
1052
1053 if (hasError)
1054 {
1055 if (debug || report)
1056 {
1057 Info<< " ***Faces not in upper triangular order." << endl;
1058 }
1059
1060 return true;
1061 }
1062
1063 if (debug || report)
1064 {
1065 Info<< " Upper triangular ordering OK." << endl;
1066 }
1067
1068 return false;
1069}
1070
1071
1073(
1074 const bool report,
1075 labelHashSet* setPtr
1076) const
1077{
1078 DebugInFunction << "Checking topological cell openness" << endl;
1079
1080 label nOpenCells = 0;
1081
1082 const faceList& f = faces();
1083 const cellList& c = cells();
1084
1085 forAll(c, celli)
1086 {
1087 const labelList& curFaces = c[celli];
1088
1089 const edgeList cellEdges = c[celli].edges(f);
1090
1091 labelList edgeUsage(cellEdges.size(), Zero);
1092
1093 forAll(curFaces, facei)
1094 {
1095 edgeList curFaceEdges = f[curFaces[facei]].edges();
1096
1097 forAll(curFaceEdges, faceEdgeI)
1098 {
1099 const edge& curEdge = curFaceEdges[faceEdgeI];
1100
1101 forAll(cellEdges, cellEdgeI)
1102 {
1103 if (cellEdges[cellEdgeI] == curEdge)
1104 {
1105 edgeUsage[cellEdgeI]++;
1106 break;
1107 }
1108 }
1109 }
1110 }
1111
1112 edgeList singleEdges(cellEdges.size());
1113 label nSingleEdges = 0;
1114
1115 forAll(edgeUsage, edgeI)
1116 {
1117 if (edgeUsage[edgeI] == 1)
1118 {
1119 singleEdges[nSingleEdges] = cellEdges[edgeI];
1120 nSingleEdges++;
1121 }
1122 else if (edgeUsage[edgeI] != 2)
1123 {
1124 if (setPtr)
1125 {
1126 setPtr->insert(celli);
1127 }
1128 }
1129 }
1130
1131 if (nSingleEdges > 0)
1132 {
1133 if (setPtr)
1134 {
1135 setPtr->insert(celli);
1136 }
1137
1138 nOpenCells++;
1139 }
1140 }
1141
1142 reduce(nOpenCells, sumOp<label>());
1143
1144 if (nOpenCells > 0)
1145 {
1146 if (debug || report)
1147 {
1148 Info<< " ***Open cells found, number of cells: " << nOpenCells
1149 << ". This problem may be fixable using the zipUpMesh utility."
1150 << endl;
1151 }
1152
1153 return true;
1154 }
1155
1156 if (debug || report)
1157 {
1158 Info<< " Topological cell zip-up check OK." << endl;
1159 }
1160
1161 return false;
1162}
1163
1164
1166(
1167 const bool report,
1168 labelHashSet* setPtr
1169) const
1170{
1171 DebugInFunction << "Checking face vertices" << endl;
1172
1173 // Check that all vertex labels are valid
1174 const faceList& f = faces();
1175
1176 label nErrorFaces = 0;
1177
1178 forAll(f, fI)
1179 {
1180 const face& curFace = f[fI];
1181
1182 if (min(curFace) < 0 || max(curFace) > nPoints())
1183 {
1184 if (setPtr)
1185 {
1186 setPtr->insert(fI);
1187 }
1188
1189 nErrorFaces++;
1190 }
1191
1192 // Uniqueness of vertices
1193 labelHashSet facePoints(2*curFace.size());
1194
1195 forAll(curFace, fp)
1196 {
1197 bool inserted = facePoints.insert(curFace[fp]);
1198
1199 if (!inserted)
1200 {
1201 if (setPtr)
1202 {
1203 setPtr->insert(fI);
1204 }
1205
1206 nErrorFaces++;
1207 }
1208 }
1209 }
1210
1211 reduce(nErrorFaces, sumOp<label>());
1212
1213 if (nErrorFaces > 0)
1214 {
1215 if (debug || report)
1216 {
1217 Info<< " Faces with invalid vertex labels found, "
1218 << " number of faces: " << nErrorFaces << endl;
1219 }
1220
1221 return true;
1222 }
1223
1224 if (debug || report)
1225 {
1226 Info<< " Face vertices OK." << endl;
1227 }
1228
1229 return false;
1230}
1231
1232
1234(
1235 const bool report,
1236 labelHashSet* setPtr
1237) const
1238{
1239 DebugInFunction << "Checking points" << endl;
1240
1241 label nFaceErrors = 0;
1242 label nCellErrors = 0;
1243
1244 const labelListList& pf = pointFaces();
1245
1246 forAll(pf, pointi)
1247 {
1248 if (pf[pointi].empty())
1249 {
1250 if (setPtr)
1251 {
1252 setPtr->insert(pointi);
1253 }
1254
1255 nFaceErrors++;
1256 }
1257 }
1258
1259
1260 forAll(pf, pointi)
1261 {
1262 const labelList& pc = pointCells(pointi);
1263
1264 if (pc.empty())
1265 {
1266 if (setPtr)
1267 {
1268 setPtr->insert(pointi);
1269 }
1270
1271 nCellErrors++;
1272 }
1273 }
1274
1275 reduce(nFaceErrors, sumOp<label>());
1276 reduce(nCellErrors, sumOp<label>());
1277
1278 if (nFaceErrors > 0 || nCellErrors > 0)
1279 {
1280 if (debug || report)
1281 {
1282 Info<< " ***Unused points found in the mesh, "
1283 "number unused by faces: " << nFaceErrors
1284 << " number unused by cells: " << nCellErrors
1285 << endl;
1286 }
1287
1288 return true;
1289 }
1290
1291 if (debug || report)
1292 {
1293 Info<< " Point usage OK." << endl;
1294 }
1295
1296 return false;
1297}
1298
1299
1301(
1302 const label facei,
1303 const Map<label>& nCommonPoints,
1304 label& nBaffleFaces,
1305 labelHashSet* setPtr
1306) const
1307{
1308 bool error = false;
1309
1310 forAllConstIters(nCommonPoints, iter)
1311 {
1312 label nbFacei = iter.key();
1313 label nCommon = iter();
1314
1315 const face& curFace = faces()[facei];
1316 const face& nbFace = faces()[nbFacei];
1317
1318 if (nCommon == nbFace.size() || nCommon == curFace.size())
1319 {
1320 if (nbFace.size() != curFace.size())
1321 {
1322 error = true;
1323 }
1324 else
1325 {
1326 nBaffleFaces++;
1327 }
1328
1329 if (setPtr)
1330 {
1331 setPtr->insert(facei);
1332 setPtr->insert(nbFacei);
1333 }
1335 }
1336
1337 return error;
1338}
1339
1340
1342(
1343 const label facei,
1344 const Map<label>& nCommonPoints,
1345 labelHashSet* setPtr
1346) const
1347{
1348 bool error = false;
1349
1350 forAllConstIters(nCommonPoints, iter)
1351 {
1352 label nbFacei = iter.key();
1353 label nCommon = iter();
1354
1355 const face& curFace = faces()[facei];
1356 const face& nbFace = faces()[nbFacei];
1357
1358 if
1359 (
1360 nCommon >= 2
1361 && nCommon != nbFace.size()
1362 && nCommon != curFace.size()
1363 )
1364 {
1365 forAll(curFace, fp)
1366 {
1367 // Get the index in the neighbouring face shared with curFace
1368 label nb = nbFace.find(curFace[fp]);
1369
1370 if (nb != -1)
1371 {
1372
1373 // Check the whole face from nb onwards for shared vertices
1374 // with neighbouring face. Rule is that any shared vertices
1375 // should be consecutive on both faces i.e. if they are
1376 // vertices fp,fp+1,fp+2 on one face they should be
1377 // vertices nb, nb+1, nb+2 (or nb+2, nb+1, nb) on the
1378 // other face.
1379
1380
1381 // Vertices before and after on curFace
1382 label fpPlus1 = curFace.fcIndex(fp);
1383 label fpMin1 = curFace.rcIndex(fp);
1384
1385 // Vertices before and after on nbFace
1386 label nbPlus1 = nbFace.fcIndex(nb);
1387 label nbMin1 = nbFace.rcIndex(nb);
1388
1389 // Find order of walking by comparing next points on both
1390 // faces.
1391 label curInc = labelMax;
1392 label nbInc = labelMax;
1393
1394 if (nbFace[nbPlus1] == curFace[fpPlus1])
1395 {
1396 curInc = 1;
1397 nbInc = 1;
1398 }
1399 else if (nbFace[nbPlus1] == curFace[fpMin1])
1400 {
1401 curInc = -1;
1402 nbInc = 1;
1403 }
1404 else if (nbFace[nbMin1] == curFace[fpMin1])
1405 {
1406 curInc = -1;
1407 nbInc = -1;
1408 }
1409 else
1410 {
1411 curInc = 1;
1412 nbInc = -1;
1413 }
1414
1415
1416 // Pass1: loop until start of common vertices found.
1417 label curNb = nb;
1418 label curFp = fp;
1419
1420 do
1421 {
1422 curFp += curInc;
1423
1424 if (curFp >= curFace.size())
1425 {
1426 curFp = 0;
1427 }
1428 else if (curFp < 0)
1429 {
1430 curFp = curFace.size()-1;
1431 }
1432
1433 curNb += nbInc;
1434
1435 if (curNb >= nbFace.size())
1436 {
1437 curNb = 0;
1438 }
1439 else if (curNb < 0)
1440 {
1441 curNb = nbFace.size()-1;
1442 }
1443 } while (curFace[curFp] == nbFace[curNb]);
1444
1445
1446 // Pass2: check equality walking from curFp, curNb
1447 // in opposite order.
1448
1449 curInc = -curInc;
1450 nbInc = -nbInc;
1451
1452 for (label commonI = 0; commonI < nCommon; commonI++)
1453 {
1454 curFp += curInc;
1455
1456 if (curFp >= curFace.size())
1457 {
1458 curFp = 0;
1459 }
1460 else if (curFp < 0)
1461 {
1462 curFp = curFace.size()-1;
1463 }
1464
1465 curNb += nbInc;
1466
1467 if (curNb >= nbFace.size())
1468 {
1469 curNb = 0;
1470 }
1471 else if (curNb < 0)
1472 {
1473 curNb = nbFace.size()-1;
1474 }
1475
1476 if (curFace[curFp] != nbFace[curNb])
1477 {
1478 if (setPtr)
1479 {
1480 setPtr->insert(facei);
1481 setPtr->insert(nbFacei);
1482 }
1483
1484 error = true;
1485
1486 break;
1487 }
1488 }
1489
1490
1491 // Done the curFace - nbFace combination.
1492 break;
1493 }
1494 }
1496 }
1497
1498 return error;
1499}
1500
1501
1503(
1504 const bool report,
1505 labelHashSet* setPtr
1506) const
1507{
1508 DebugInFunction << "Checking face-face connectivity" << endl;
1509
1510 const labelListList& pf = pointFaces();
1511
1512 label nBaffleFaces = 0;
1513 label nErrorDuplicate = 0;
1514 label nErrorOrder = 0;
1515 Map<label> nCommonPoints;
1516
1517 for (label facei = 0; facei < nFaces(); facei++)
1518 {
1519 const face& curFace = faces()[facei];
1520
1521 // Calculate number of common points between current facei and
1522 // neighbouring face. Store on map.
1523 nCommonPoints.clear();
1524
1525 forAll(curFace, fp)
1526 {
1527 label pointi = curFace[fp];
1528
1529 const labelList& nbs = pf[pointi];
1530
1531 forAll(nbs, nbI)
1532 {
1533 label nbFacei = nbs[nbI];
1534
1535 if (facei < nbFacei)
1536 {
1537 // Only check once for each combination of two faces.
1538 ++(nCommonPoints(nbFacei, 0));
1539 }
1540 }
1541 }
1542
1543 // Perform various checks on common points
1544
1545 // Check all vertices shared (duplicate point)
1546 if (checkDuplicateFaces(facei, nCommonPoints, nBaffleFaces, setPtr))
1547 {
1548 nErrorDuplicate++;
1549 }
1550
1551 // Check common vertices are consecutive on both faces
1552 if (checkCommonOrder(facei, nCommonPoints, setPtr))
1553 {
1554 nErrorOrder++;
1555 }
1556 }
1557
1558 reduce(nBaffleFaces, sumOp<label>());
1559 reduce(nErrorDuplicate, sumOp<label>());
1560 reduce(nErrorOrder, sumOp<label>());
1561
1562 if (nBaffleFaces)
1563 {
1564 Info<< " Number of identical duplicate faces (baffle faces): "
1565 << nBaffleFaces << endl;
1566 }
1567
1568 if (nErrorDuplicate > 0 || nErrorOrder > 0)
1569 {
1570 // These are actually warnings, not errors.
1571 if (nErrorDuplicate > 0)
1572 {
1573 Info<< " <<Number of duplicate (not baffle) faces found: "
1574 << nErrorDuplicate
1575 << ". This might indicate a problem." << endl;
1576 }
1577
1578 if (nErrorOrder > 0)
1579 {
1580 Info<< " <<Number of faces with non-consecutive shared points: "
1581 << nErrorOrder << ". This might indicate a problem." << endl;
1582 }
1583
1584 return false; //return true;
1585 }
1586
1587 if (debug || report)
1588 {
1589 Info<< " Face-face connectivity OK." << endl;
1590 }
1592 return false;
1593}
1594
1595
1596// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1598bool Foam::primitiveMesh::checkClosedBoundary(const bool report) const
1599{
1600 return checkClosedBoundary(faceAreas(), report, bitSet());
1601}
1602
1603
1605(
1606 const bool report,
1607 labelHashSet* setPtr,
1608 labelHashSet* aspectSetPtr,
1609 const Vector<label>& solutionD
1610) const
1611{
1612 return checkClosedCells
1613 (
1614 faceAreas(),
1615 cellVolumes(),
1616 report,
1617 setPtr,
1618 aspectSetPtr,
1619 solutionD
1620 );
1621}
1622
1623
1625(
1626 const bool report,
1627 labelHashSet* setPtr
1628) const
1629{
1630 return checkFaceAreas
1631 (
1632 faceAreas(),
1633 report,
1634 false, // detailedReport,
1635 setPtr
1636 );
1637}
1638
1639
1641(
1642 const bool report,
1643 labelHashSet* setPtr
1644) const
1645{
1646 return checkCellVolumes
1647 (
1648 cellVolumes(),
1649 report,
1650 false, // detailedReport,
1651 setPtr
1652 );
1653}
1654
1655
1657(
1658 const bool report,
1659 labelHashSet* setPtr
1660) const
1661{
1662 return checkFaceOrthogonality
1663 (
1664 faceAreas(),
1666 report,
1667 setPtr
1668 );
1669}
1670
1671
1673(
1674 const bool report,
1675 const scalar minPyrVol,
1676 labelHashSet* setPtr
1677) const
1678{
1679 return checkFacePyramids
1680 (
1681 points(),
1682 cellCentres(),
1683 report,
1684 false, // detailedReport,
1685 minPyrVol,
1686 setPtr
1687 );
1688}
1689
1690
1692(
1693 const bool report,
1694 labelHashSet* setPtr
1695) const
1696{
1697 return checkFaceSkewness
1698 (
1699 points(),
1700 faceCentres(),
1701 faceAreas(),
1703 report,
1704 setPtr
1705 );
1706}
1707
1708
1710(
1711 const bool report,
1712 const scalar maxDeg,
1713 labelHashSet* setPtr
1714) const
1715{
1716 return checkFaceAngles
1717 (
1718 points(),
1719 faceAreas(),
1720 report,
1721 maxDeg,
1722 setPtr
1723 );
1724}
1725
1726
1728(
1729 const bool report,
1730 const scalar warnFlatness,
1731 labelHashSet* setPtr
1732) const
1733{
1734 return checkFaceFlatness
1735 (
1736 points(),
1737 faceCentres(),
1738 faceAreas(),
1739 report,
1740 warnFlatness,
1741 setPtr
1742 );
1743}
1744
1745
1747(
1748 const bool report,
1749 labelHashSet* setPtr
1750) const
1751{
1752 return checkConcaveCells
1753 (
1754 faceAreas(),
1756 report,
1757 setPtr
1758 );
1759}
1760
1761
1762bool Foam::primitiveMesh::checkTopology(const bool report) const
1763{
1764 label nFailedChecks = 0;
1765
1766 if (checkPoints(report)) ++nFailedChecks;
1767 if (checkUpperTriangular(report)) ++nFailedChecks;
1768 if (checkCellsZipUp(report)) ++nFailedChecks;
1769 if (checkFaceVertices(report)) ++nFailedChecks;
1770 if (checkFaceFaces(report)) ++nFailedChecks;
1771
1772 if (nFailedChecks)
1773 {
1774 if (debug || report)
1775 {
1776 Info<< " Failed " << nFailedChecks
1777 << " mesh topology checks." << endl;
1778 }
1779
1780 return true;
1781 }
1782
1783 if (debug || report)
1784 {
1785 Info<< " Mesh topology OK." << endl;
1786 }
1787
1788 return false;
1789}
1790
1791
1792bool Foam::primitiveMesh::checkGeometry(const bool report) const
1793{
1794 label nFailedChecks = 0;
1795
1796 if (checkClosedBoundary(report)) ++nFailedChecks;
1797 if (checkClosedCells(report)) ++nFailedChecks;
1798 if (checkFaceAreas(report)) ++nFailedChecks;
1799 if (checkCellVolumes(report)) ++nFailedChecks;
1800 if (checkFaceOrthogonality(report)) ++nFailedChecks;
1801 if (checkFacePyramids(report)) ++nFailedChecks;
1802 if (checkFaceSkewness(report)) ++nFailedChecks;
1803
1804 if (nFailedChecks)
1805 {
1806 if (debug || report)
1807 {
1808 Info<< " Failed " << nFailedChecks
1809 << " mesh geometry checks." << endl;
1810 }
1811
1812 return true;
1813 }
1814
1815 if (debug || report)
1816 {
1817 Info<< " Mesh geometry OK." << endl;
1818 }
1819
1820 return false;
1821}
1822
1823
1824bool Foam::primitiveMesh::checkMesh(const bool report) const
1825{
1826 DebugInFunction << "Checking primitiveMesh" << endl;
1827
1828 label nFailedChecks = checkTopology(report) + checkGeometry(report);
1829
1830 if (nFailedChecks)
1831 {
1832 if (debug || report)
1833 {
1834 Info<< " Failed " << nFailedChecks
1835 << " mesh checks." << endl;
1836 }
1837
1838 return true;
1839 }
1840
1841 if (debug || report)
1842 {
1843 Info<< "Mesh OK." << endl;
1844 }
1845
1846 return false;
1847}
1848
1849
1850Foam::scalar Foam::primitiveMesh::setClosedThreshold(const scalar val)
1851{
1852 scalar prev = closedThreshold_;
1853 closedThreshold_ = val;
1854
1855 return prev;
1856}
1857
1858
1859Foam::scalar Foam::primitiveMesh::setAspectThreshold(const scalar val)
1860{
1861 scalar prev = aspectThreshold_;
1862 aspectThreshold_ = val;
1863
1864 return prev;
1865}
1866
1867
1868Foam::scalar Foam::primitiveMesh::setNonOrthThreshold(const scalar val)
1869{
1870 scalar prev = nonOrthThreshold_;
1871 nonOrthThreshold_ = val;
1872
1873 return prev;
1874}
1875
1876
1877Foam::scalar Foam::primitiveMesh::setSkewThreshold(const scalar val)
1878{
1879 scalar prev = skewThreshold_;
1880 skewThreshold_ = val;
1881
1882 return prev;
1883}
1884
1885
1886// ************************************************************************* //
Various functions to operate on Lists.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
void clear()
Remove all entries from table.
Definition HashTable.C:742
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
label size() const noexcept
Number of entries.
Definition PackedList.H:392
A list that is sorted upon construction or when explicitly requested with the sort() method.
const labelList & indices() const noexcept
Return the list of sorted indices. Updated every sort.
void sort()
Forward (stable) sort the list (if changed after construction).
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
Definition UListI.H:106
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition UListI.H:99
static void reduceOr(bool &value, const int communicator=worldComm)
Logical (or) reduction (MPI_AllReduce).
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition error.H:74
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
label find(const Foam::edge &e) const
Find the edge within the face.
Definition face.C:791
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition pointCells.H:55
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) most concave angle between two conse...
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 tmp< scalarField > faceSkewness(const primitiveMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
static tmp< scalarField > faceOrthogonality(const primitiveMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate non-orthogonality field (internal faces only).
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' normals against the face average nor...
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
bool checkFaceOrthogonality(const vectorField &fAreas, const vectorField &cellCtrs, const bool report, labelHashSet *setPtr) const
Check for non-orthogonality.
static scalar planarCosAngle_
Threshold where faces are considered coplanar.
virtual const labelList & faceOwner() const =0
Face face-owner addressing.
const labelListList & cellEdges() const
virtual const labelList & faceNeighbour() const =0
Face face-neighbour addressing.
const vectorField & faceCentres() const
virtual bool checkCellsZipUp(const bool report=false, labelHashSet *setPtr=nullptr) const
Check cell zip-up.
static scalar skewThreshold_
Skewness warning threshold.
static scalar setNonOrthThreshold(const scalar)
Set the non-orthogonality warning threshold in degrees.
const labelListList & pointCells() const
const scalarField & cellVolumes() const
virtual const faceList & faces() const =0
Return faces.
label nInternalFaces() const noexcept
Number of internal faces.
virtual bool checkMesh(const bool report=false) const
Check mesh for correctness. Returns false for no error.
virtual bool checkTopology(const bool report=false) const
Check mesh topology for correctness.
bool checkCommonOrder(const label, const Map< label > &, labelHashSet *) const
Check that shared points are in consecutive order.
virtual bool checkUpperTriangular(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face ordering.
static scalar closedThreshold_
Static data to control mesh checking.
const vectorField & cellCentres() const
bool checkFacePyramids(const pointField &points, const vectorField &ctrs, const bool report, const bool detailedReport, const scalar minPyrVol, labelHashSet *setPtr) const
Check face pyramid volume.
bool checkFaceAngles(const pointField &points, const vectorField &faceAreas, const bool report, const scalar maxDeg, labelHashSet *setPtr) const
Check face angles.
static scalar nonOrthThreshold_
Non-orthogonality warning threshold in deg.
static scalar setSkewThreshold(const scalar)
Set the skewness warning threshold as percentage.
static scalar setAspectThreshold(const scalar)
Set the aspect ratio warning threshold.
bool checkDuplicateFaces(const label, const Map< label > &, label &nBaffleFaces, labelHashSet *) const
Check if all points on face are shared with another face.
virtual bool checkFaceVertices(const bool report=false, labelHashSet *setPtr=nullptr) const
Check uniqueness of face vertices.
bool checkCellVolumes(const scalarField &vols, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative cell volumes.
bool checkFaceAreas(const vectorField &faceAreas, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative face areas.
label nFaces() const noexcept
Number of mesh faces.
virtual bool checkPoints(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for unused points.
bool checkFaceSkewness(const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs, const bool report, labelHashSet *setPtr) const
Check face skewness.
bool checkFaceFlatness(const pointField &points, const vectorField &faceCentres, const vectorField &faceAreas, const bool report, const scalar warnFlatness, labelHashSet *setPtr) const
Check face warpage.
const labelListList & pointFaces() const
static scalar aspectThreshold_
Aspect ratio warning threshold.
const vectorField & faceAreas() const
bool checkConcaveCells(const vectorField &fAreas, const pointField &fCentres, const bool report, labelHashSet *setPtr) const
Check for concave cells by the planes of faces.
bool checkClosedBoundary(const vectorField &areas, const bool report, const bitSet &internalOrCoupledFaces) const
Check boundary for closedness.
virtual bool checkFaceFaces(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face-face connectivity.
virtual const pointField & points() const =0
Return mesh points.
bool checkClosedCells(const vectorField &faceAreas, const scalarField &cellVolumes, const bool report, labelHashSet *setPtr, labelHashSet *aspectSetPtr, const Vector< label > &meshD) const
Check cells for closedness.
virtual bool checkGeometry(const bool report=false) const
Check mesh geometry (& implicitly topology) for correctness.
static scalar setClosedThreshold(const scalar)
Set the closedness ratio warning threshold.
A class for managing temporary objects.
Definition tmp.H:75
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
label nPoints
const cellShapeList & cells
#define DebugInFunction
Report an information message using Foam::Info.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
List< edge > edgeList
List of edge.
Definition edgeList.H:32
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar asin(const dimensionedScalar &ds)
Type gSum(const FieldField< Field, Type > &f)
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
label checkGeometry(const polyMesh &mesh, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, autoPtr< coordSetWriter > &setWriter)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
label checkTopology(const polyMesh &mesh, const bool allTopology, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, autoPtr< coordSetWriter > &setWriter, const bool writeBadEdges=false)
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Vector< solveScalar > solveVector
Definition vector.H:60
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Unit conversion functions.