polyMeshGeometry.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-2017 OpenFOAM Foundation
9 Copyright (C) 2015-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 "polyMeshGeometry.H"
31#include "pyramidPointFaceRef.H"
32#include "tetPointRef.H"
33#include "syncTools.H"
34#include "unitConversion.H"
35#include "primitiveMeshTools.H"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace Foam
40{
42}
43
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
47void Foam::polyMeshGeometry::updateFaceCentresAndAreas
48(
49 const pointField& p,
50 const labelList& changedFaces
51)
52{
53 const faceList& fcs = mesh_.faces();
54
55 for (const label facei : changedFaces)
56 {
57 const face& f = fcs[facei];
58 const label nPoints = f.size();
59
60 // If the face is a triangle, do a direct calculation for efficiency
61 // and to avoid round-off error-related problems
62 if (nPoints == 3)
63 {
64 faceCentres_[facei] = (1.0/3.0)*(p[f[0]] + p[f[1]] + p[f[2]]);
65 faceAreas_[facei] = 0.5*((p[f[1]] - p[f[0]])^(p[f[2]] - p[f[0]]));
66 }
67 else
68 {
69 solveVector sumN = Zero;
70 solveScalar sumA = Zero;
71 solveVector sumAc = Zero;
72
73 solveVector fCentre = p[f[0]];
74 for (label pi = 1; pi < nPoints; ++pi)
75 {
76 fCentre += solveVector(p[f[pi]]);
77 }
78 fCentre /= nPoints;
79
80 for (label pi = 0; pi < nPoints; ++pi)
81 {
82 const solveVector thisPoint(p[f.thisLabel(pi)]);
83 const solveVector nextPoint(p[f.nextLabel(pi)]);
84
85 solveVector c = thisPoint + nextPoint + fCentre;
86 solveVector n = (nextPoint - thisPoint)^(fCentre - thisPoint);
87 solveScalar a = mag(n);
88
89 sumN += n;
90 sumA += a;
91 sumAc += a*c;
92 }
93
94 if (sumA < ROOTVSMALL)
95 {
96 faceCentres_[facei] = fCentre;
97 faceAreas_[facei] = Zero;
98 }
99 else
100 {
101 faceCentres_[facei] = (1.0/3.0)*sumAc/sumA;
102 faceAreas_[facei] = 0.5*sumN;
103 }
104 }
105 }
106}
107
108
109void Foam::polyMeshGeometry::updateCellCentresAndVols
110(
111 const labelList& changedCells,
112 const labelList& changedFaces
113)
114{
115 const labelList& own = mesh().faceOwner();
116 const cellList& cells = mesh().cells();
117
118 // Clear the fields for accumulation
119 UIndirectList<vector>(cellCentres_, changedCells) = Zero;
120 UIndirectList<scalar>(cellVolumes_, changedCells) = Zero;
121
122
123 // Re-calculate the changed cell centres and volumes
124 for (const label celli : changedCells)
125 {
126 const labelList& cFaces = cells[celli];
127
128 // Estimate the cell centre and bounding box using the face centres
129 vector cEst(Zero);
130 boundBox bb(boundBox::invertedBox);
131
132 for (const label facei : cFaces)
133 {
134 const point& fc = faceCentres_[facei];
135 cEst += fc;
136 bb.add(fc);
137 }
138 cEst /= cFaces.size();
139
140
141 // Sum up the face-pyramid contributions
142 for (const label facei : cFaces)
143 {
144 // Calculate 3* the face-pyramid volume
145 scalar pyr3Vol = faceAreas_[facei] & (faceCentres_[facei] - cEst);
146
147 if (own[facei] != celli)
148 {
149 pyr3Vol = -pyr3Vol;
150 }
151
152 // Accumulate face-pyramid volume
153 cellVolumes_[celli] += pyr3Vol;
154
155 // Calculate the face-pyramid centre
156 const vector pCtr = (3.0/4.0)*faceCentres_[facei] + (1.0/4.0)*cEst;
157
158 // Accumulate volume-weighted face-pyramid centre
159 cellCentres_[celli] += pyr3Vol*pCtr;
160 }
161
162 // Average the accumulated quantities
163
164 if (mag(cellVolumes_[celli]) > VSMALL)
165 {
166 point cc = cellCentres_[celli] / cellVolumes_[celli];
167
168 // Do additional check for collapsed cells since some volumes
169 // (e.g. 1e-33) do not trigger above but do return completely
170 // wrong cell centre
171 if (bb.contains(cc))
172 {
173 cellCentres_[celli] = cc;
174 }
175 else
176 {
177 cellCentres_[celli] = cEst;
178 }
179 }
180 else
181 {
182 cellCentres_[celli] = cEst;
183 }
184
185 cellVolumes_[celli] *= (1.0/3.0);
186 }
187}
188
189
191(
192 const polyMesh& mesh,
193 const labelList& changedFaces
194)
195{
196 const labelList& own = mesh.faceOwner();
197 const labelList& nei = mesh.faceNeighbour();
198
199 labelHashSet affectedCells(2*changedFaces.size());
200
201 for (const label facei : changedFaces)
202 {
203 affectedCells.insert(own[facei]);
204
205 if (mesh.isInternalFace(facei))
206 {
207 affectedCells.insert(nei[facei]);
208 }
209 }
210 return affectedCells.toc();
211}
212
213
214Foam::scalar Foam::polyMeshGeometry::checkNonOrtho
215(
216 const polyMesh& mesh,
217 const bool report,
218 const scalar severeNonorthogonalityThreshold,
219 const label facei,
220 const vector& s, // face area vector
221 const vector& d, // cc-cc vector
222
223 label& severeNonOrth,
224 label& errorNonOrth,
225 labelHashSet* setPtr
226)
227{
228 scalar dDotS = (d & s)/(mag(d)*mag(s) + VSMALL);
229
230 if (dDotS < severeNonorthogonalityThreshold)
231 {
232 label nei = -1;
233
234 if (mesh.isInternalFace(facei))
235 {
236 nei = mesh.faceNeighbour()[facei];
237 }
238
239 if (dDotS > SMALL)
240 {
241 if (report)
242 {
243 // Severe non-orthogonality but mesh still OK
244 Pout<< "Severe non-orthogonality for face " << facei
245 << " between cells " << mesh.faceOwner()[facei]
246 << " and " << nei
247 << ": Angle = "
248 << radToDeg(::acos(dDotS))
249 << " deg." << endl;
250 }
251
252 ++severeNonOrth;
253 }
254 else
255 {
256 // Non-orthogonality greater than 90 deg
257 if (report)
258 {
260 << "Severe non-orthogonality detected for face "
261 << facei
262 << " between cells " << mesh.faceOwner()[facei]
263 << " and " << nei
264 << ": Angle = "
265 << radToDeg(::acos(dDotS))
266 << " deg." << endl;
267 }
268
269 ++errorNonOrth;
270 }
271
272 if (setPtr)
273 {
274 setPtr->insert(facei);
275 }
276 }
277 return dDotS;
278}
279
280
281// Create the neighbour pyramid - it will have positive volume
282bool Foam::polyMeshGeometry::checkFaceTet
283(
284 const polyMesh& mesh,
285 const bool report,
286 const scalar minTetQuality,
287 const pointField& p,
288 const label facei,
289 const point& fc, // face centre
290 const point& cc, // cell centre
291
292 labelHashSet* setPtr
293)
294{
295 const face& f = mesh.faces()[facei];
296
297 forAll(f, fp)
298 {
299 scalar tetQual = tetPointRef
300 (
301 p[f[fp]],
302 p[f.nextLabel(fp)],
303 fc,
304 cc
305 ).quality();
306
307 if (tetQual < minTetQuality)
308 {
309 if (report)
310 {
311 Pout<< "bool polyMeshGeometry::checkFaceTets("
312 << "const bool, const scalar, const pointField&"
313 << ", const pointField&"
314 << ", const labelList&, labelHashSet*) : "
315 << "face " << facei
316 << " has a triangle that points the wrong way." << nl
317 << "Tet quality: " << tetQual
318 << " Face " << facei
319 << endl;
320 }
321 if (setPtr)
322 {
323 setPtr->insert(facei);
324 }
325 return true;
326 }
327 }
328
329 return false;
330}
331
332
333// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
334
336:
337 mesh_(mesh)
338{
339 correct();
340}
341
342
343// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
344
346{
347 faceAreas_ = mesh_.faceAreas();
348 faceCentres_ = mesh_.faceCentres();
349 cellCentres_ = mesh_.cellCentres();
350 cellVolumes_ = mesh_.cellVolumes();
351}
352
353
355(
356 const pointField& p,
357 const labelList& changedFaces
358)
359{
360 // Update face quantities
361 updateFaceCentresAndAreas(p, changedFaces);
362 // Update cell quantities from face quantities
363 updateCellCentresAndVols(affectedCells(mesh_, changedFaces), changedFaces);
364}
365
366
368(
369 const bool report,
370 const scalar orthWarn,
371 const polyMesh& mesh,
372 const vectorField& cellCentres,
373 const vectorField& faceAreas,
374 const labelList& checkFaces,
375 const List<labelPair>& baffles,
376 labelHashSet* setPtr
377)
378{
379 // for all internal and coupled faces check that the d dot S product
380 // is positive
381
382 const labelList& own = mesh.faceOwner();
383 const labelList& nei = mesh.faceNeighbour();
385
386 // Severe nonorthogonality threshold
387 const scalar severeNonorthogonalityThreshold = ::cos(degToRad(orthWarn));
388
389 // Calculate coupled cell centre
391
392 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
393 {
394 neiCc[facei-mesh.nInternalFaces()] = cellCentres[own[facei]];
395 }
396
398
399 scalar minDDotS = GREAT;
400
401 scalar sumDDotS = 0;
402 label nDDotS = 0;
403
404 label severeNonOrth = 0;
405
406 label errorNonOrth = 0;
407
408 for (const label facei : checkFaces)
409 {
410 const point& ownCc = cellCentres[own[facei]];
411
412 if (mesh.isInternalFace(facei))
413 {
414 scalar dDotS = checkNonOrtho
415 (
416 mesh,
417 report,
418 severeNonorthogonalityThreshold,
419 facei,
420 faceAreas[facei],
421 cellCentres[nei[facei]] - ownCc,
422
423 severeNonOrth,
424 errorNonOrth,
425 setPtr
426 );
427
428 if (dDotS < minDDotS)
429 {
430 minDDotS = dDotS;
431 }
432
433 sumDDotS += dDotS;
434 ++nDDotS;
435 }
436 else
437 {
438 const label patchi = patches.whichPatch(facei);
439
440 if (patches[patchi].coupled())
441 {
442 scalar dDotS = checkNonOrtho
443 (
444 mesh,
445 report,
446 severeNonorthogonalityThreshold,
447 facei,
448 faceAreas[facei],
449 neiCc[facei-mesh.nInternalFaces()] - ownCc,
450
451 severeNonOrth,
452 errorNonOrth,
453 setPtr
454 );
455
456 if (dDotS < minDDotS)
457 {
458 minDDotS = dDotS;
459 }
460
461 sumDDotS += dDotS;
462 ++nDDotS;
463 }
464 }
465 }
466
467 for (const labelPair& baffle : baffles)
468 {
469 const label face0 = baffle.first();
470 const label face1 = baffle.second();
471
472 const point& ownCc = cellCentres[own[face0]];
473
474 scalar dDotS = checkNonOrtho
475 (
476 mesh,
477 report,
478 severeNonorthogonalityThreshold,
479 face0,
480 faceAreas[face0],
481 cellCentres[own[face1]] - ownCc,
482
483 severeNonOrth,
484 errorNonOrth,
485 setPtr
486 );
487
488 if (dDotS < minDDotS)
489 {
490 minDDotS = dDotS;
491 }
492
493 sumDDotS += dDotS;
494 ++nDDotS;
495 }
496
497 reduce(minDDotS, minOp<scalar>());
498 reduce(sumDDotS, sumOp<scalar>());
499 reduce(nDDotS, sumOp<label>());
500 reduce(severeNonOrth, sumOp<label>());
501 reduce(errorNonOrth, sumOp<label>());
502
503 // Only report if there are some internal faces
504 if (nDDotS > 0)
505 {
506 if (report && minDDotS < severeNonorthogonalityThreshold)
507 {
508 Info<< "Number of non-orthogonality errors: " << errorNonOrth
509 << ". Number of severely non-orthogonal faces: "
510 << severeNonOrth << "." << endl;
511 }
512 }
513
514 if (report)
515 {
516 if (nDDotS > 0)
517 {
518 Info<< "Mesh non-orthogonality Max: "
519 << radToDeg(::acos(minDDotS))
520 << " average: " << radToDeg(::acos(sumDDotS/nDDotS))
521 << endl;
522 }
523 }
524
525 if (errorNonOrth > 0)
526 {
527 if (report)
528 {
530 << "Error in non-orthogonality detected" << endl;
531 }
532
533 return true;
534 }
535
536 if (report)
537 {
538 Info<< "Non-orthogonality check OK.\n" << endl;
539 }
540
541 return false;
542}
543
544
546(
547 const bool report,
548 const scalar minPyrVol,
549 const polyMesh& mesh,
550 const vectorField& cellCentres,
551 const pointField& p,
552 const labelList& checkFaces,
553 const List<labelPair>& baffles,
554 labelHashSet* setPtr
555)
556{
557 // check whether face area vector points to the cell with higher label
558 const labelList& own = mesh.faceOwner();
559 const labelList& nei = mesh.faceNeighbour();
560
561 const faceList& f = mesh.faces();
562
563 label nErrorPyrs = 0;
564
565 for (const label facei : checkFaces)
566 {
567 // Create the owner pyramid - it will have negative volume
568 scalar pyrVol = pyramidPointFaceRef
569 (
570 f[facei],
571 cellCentres[own[facei]]
572 ).mag(p);
573
574 if (pyrVol > -minPyrVol)
575 {
576 ++nErrorPyrs;
577
578 if (report)
579 {
580 Pout<< "bool polyMeshGeometry::checkFacePyramids("
581 << "const bool, const scalar, const pointField&"
582 << ", const labelList&, labelHashSet*): "
583 << "face " << facei << " points the wrong way. " << endl
584 << "Pyramid volume: " << -pyrVol
585 << " Face " << f[facei] << " area: " << f[facei].mag(p)
586 << " Owner cell: " << own[facei] << endl
587 << "Owner cell vertex labels: "
588 << mesh.cells()[own[facei]].labels(f)
589 << endl;
590 }
591 if (setPtr)
592 {
593 setPtr->insert(facei);
594 }
595 }
596
597 if (mesh.isInternalFace(facei))
598 {
599 // Create the neighbour pyramid - it will have positive volume
600 scalar pyrVol =
601 pyramidPointFaceRef(f[facei], cellCentres[nei[facei]]).mag(p);
602
603 if (pyrVol < minPyrVol)
604 {
605 ++nErrorPyrs;
606
607 if (report)
608 {
609 Pout<< "bool polyMeshGeometry::checkFacePyramids("
610 << "const bool, const scalar, const pointField&"
611 << ", const labelList&, labelHashSet*): "
612 << "face " << facei << " points the wrong way. " << endl
613 << "Pyramid volume: " << -pyrVol
614 << " Face " << f[facei] << " area: " << f[facei].mag(p)
615 << " Neighbour cell: " << nei[facei] << endl
616 << "Neighbour cell vertex labels: "
617 << mesh.cells()[nei[facei]].labels(f)
618 << endl;
619 }
620 if (setPtr)
621 {
622 setPtr->insert(facei);
623 }
624 }
625 }
626 }
627
628 for (const labelPair& baffle : baffles)
629 {
630 const label face0 = baffle.first();
631 const label face1 = baffle.second();
632
633 const point& ownCc = cellCentres[own[face0]];
634
635 // Create the owner pyramid - it will have negative volume
636 scalar pyrVolOwn = pyramidPointFaceRef
637 (
638 f[face0],
639 ownCc
640 ).mag(p);
641
642 if (pyrVolOwn > -minPyrVol)
643 {
644 ++nErrorPyrs;
645
646 if (report)
647 {
648 Pout<< "bool polyMeshGeometry::checkFacePyramids("
649 << "const bool, const scalar, const pointField&"
650 << ", const labelList&, labelHashSet*): "
651 << "face " << face0 << " points the wrong way. " << endl
652 << "Pyramid volume: " << -pyrVolOwn
653 << " Face " << f[face0] << " area: " << f[face0].mag(p)
654 << " Owner cell: " << own[face0] << endl
655 << "Owner cell vertex labels: "
656 << mesh.cells()[own[face0]].labels(f)
657 << endl;
658 }
659 if (setPtr)
660 {
661 setPtr->insert(face0);
662 }
663 }
664
665 // Create the neighbour pyramid - it will have positive volume
666 scalar pyrVolNbr =
667 pyramidPointFaceRef(f[face0], cellCentres[own[face1]]).mag(p);
668
669 if (pyrVolNbr < minPyrVol)
670 {
671 ++nErrorPyrs;
672
673 if (report)
674 {
675 Pout<< "bool polyMeshGeometry::checkFacePyramids("
676 << "const bool, const scalar, const pointField&"
677 << ", const labelList&, labelHashSet*): "
678 << "face " << face0 << " points the wrong way. " << endl
679 << "Pyramid volume: " << -pyrVolNbr
680 << " Face " << f[face0] << " area: " << f[face0].mag(p)
681 << " Neighbour cell: " << own[face1] << endl
682 << "Neighbour cell vertex labels: "
683 << mesh.cells()[own[face1]].labels(f)
684 << endl;
685 }
686 if (setPtr)
687 {
688 setPtr->insert(face0);
689 }
690 }
691 }
692
693 reduce(nErrorPyrs, sumOp<label>());
694
695 if (nErrorPyrs > 0)
696 {
697 if (report)
698 {
700 << "Error in face pyramids: faces pointing the wrong way."
701 << endl;
702 }
703
704 return true;
705 }
706
707 if (report)
708 {
709 Info<< "Face pyramids OK.\n" << endl;
710 }
711
712 return false;
713}
714
715
717(
718 const bool report,
719 const scalar minTetQuality,
720 const polyMesh& mesh,
721 const vectorField& cellCentres,
722 const vectorField& faceCentres,
723 const pointField& p,
724 const labelList& checkFaces,
725 const List<labelPair>& baffles,
726 labelHashSet* setPtr
727)
728{
729 // check whether decomposing each cell into tets results in
730 // positive volume, non-flat tets
731 const labelList& own = mesh.faceOwner();
732 const labelList& nei = mesh.faceNeighbour();
734
735 // Calculate coupled cell centre
737
738 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
739 {
740 neiCc[facei - mesh.nInternalFaces()] = cellCentres[own[facei]];
741 }
742
744
745 label nErrorTets = 0;
746
747 for (const label facei : checkFaces)
748 {
749 // Create the owner pyramid - note: exchange cell and face centre
750 // to get positive volume.
751 bool tetError = checkFaceTet
752 (
753 mesh,
754 report,
755 minTetQuality,
756 p,
757 facei,
758 cellCentres[own[facei]], // face centre
759 faceCentres[facei], // cell centre
760 setPtr
761 );
762
763 if (tetError)
764 {
765 ++nErrorTets;
766 }
767
768 if (mesh.isInternalFace(facei))
769 {
770 // Create the neighbour tets - they will have positive volume
771 bool tetError = checkFaceTet
772 (
773 mesh,
774 report,
775 minTetQuality,
776 p,
777 facei,
778 faceCentres[facei], // face centre
779 cellCentres[nei[facei]], // cell centre
780 setPtr
781 );
782
783 if (tetError)
784 {
785 ++nErrorTets;
786 }
787
788 if
789 (
791 (
792 mesh,
793 facei,
794 minTetQuality,
795 report
796 ) == -1
797 )
798 {
799 ++nErrorTets;
800 if (setPtr)
801 {
802 setPtr->insert(facei);
803 }
804 }
805 }
806 else
807 {
808 label patchi = patches.whichPatch(facei);
809
810 if (patches[patchi].coupled())
811 {
812 if
813 (
815 (
816 mesh,
817 facei,
818 neiCc[facei - mesh.nInternalFaces()],
819 minTetQuality,
820 report
821 ) == -1
822 )
823 {
824 ++nErrorTets;
825 if (setPtr)
826 {
827 setPtr->insert(facei);
828 }
829 }
830 }
831 else
832 {
833 if
834 (
836 (
837 mesh,
838 facei,
839 minTetQuality,
840 report
841 ) == -1
842 )
843 {
844 ++nErrorTets;
845 if (setPtr)
846 {
847 setPtr->insert(facei);
848 }
849 }
850 }
851 }
852 }
853
854 for (const labelPair& baffle : baffles)
855 {
856 const label face0 = baffle.first();
857 const label face1 = baffle.second();
858
859 bool tetError = checkFaceTet
860 (
861 mesh,
862 report,
863 minTetQuality,
864 p,
865 face0,
866 cellCentres[own[face0]], // face centre
867 faceCentres[face0], // cell centre
868 setPtr
869 );
870
871 if (tetError)
872 {
873 ++nErrorTets;
874 }
875
876 // Create the neighbour tets - they will have positive volume
877 tetError = checkFaceTet
878 (
879 mesh,
880 report,
881 minTetQuality,
882 p,
883 face0,
884 faceCentres[face0], // face centre
885 cellCentres[own[face1]], // cell centre
886 setPtr
887 );
888
889 if (tetError)
890 {
891 ++nErrorTets;
892 }
893
894 if
895 (
897 (
898 mesh,
899 face0,
900 cellCentres[own[face1]],
901 minTetQuality,
902 report
903 ) == -1
904 )
905 {
906 ++nErrorTets;
907 if (setPtr)
908 {
909 setPtr->insert(face0);
910 }
911 }
912 }
913
914 reduce(nErrorTets, sumOp<label>());
915
916 if (nErrorTets > 0)
917 {
918 if (report)
919 {
921 << "Error in face decomposition: negative tets."
922 << endl;
923 }
924
925 return true;
926 }
927
928 if (report)
929 {
930 Info<< "Face tets OK.\n" << endl;
931 }
932
933 return false;
934}
935
936
938(
939 const bool report,
940 const scalar internalSkew,
941 const scalar boundarySkew,
942 const polyMesh& mesh,
943 const pointField& points,
944 const vectorField& cellCentres,
945 const vectorField& faceCentres,
946 const vectorField& faceAreas,
947 const labelList& checkFaces,
948 const List<labelPair>& baffles,
949 labelHashSet* setPtr
950)
951{
952 // Warn if the skew correction vector is more than skew times
953 // larger than the face area vector
954
955 const labelList& own = mesh.faceOwner();
956 const labelList& nei = mesh.faceNeighbour();
957 const faceList& faces = mesh.faces();
959
960 // Calculate coupled cell centre
961 pointField neiCc;
962 syncTools::swapBoundaryCellPositions(mesh, cellCentres, neiCc);
963
964 scalar maxSkew = 0;
965
966 label nWarnSkew = 0;
967
968 for (const label facei : checkFaces)
969 {
970 if (mesh.isInternalFace(facei))
971 {
972 scalar skewness = primitiveMeshTools::faceSkewness
973 (
974 faces,
975 points,
976 faceCentres,
977 faceAreas,
978
979 facei,
980 cellCentres[own[facei]],
981 cellCentres[nei[facei]]
982 );
983
984 // Check if the skewness vector is greater than the PN vector.
985 // This does not cause trouble but is a good indication of a poor
986 // mesh.
987 if (skewness > internalSkew)
988 {
989 ++nWarnSkew;
990
991 if (report)
992 {
993 Pout<< "Severe skewness for face " << facei
994 << " skewness = " << skewness << endl;
995 }
996 if (setPtr)
997 {
998 setPtr->insert(facei);
999 }
1000 }
1001
1002 maxSkew = max(maxSkew, skewness);
1003 }
1004 else if (patches[patches.whichPatch(facei)].coupled())
1005 {
1006 scalar skewness = primitiveMeshTools::faceSkewness
1007 (
1008 faces,
1009 points,
1010 faceCentres,
1011 faceAreas,
1012
1013 facei,
1014 cellCentres[own[facei]],
1015 neiCc[facei-mesh.nInternalFaces()]
1016 );
1017
1018 // Check if the skewness vector is greater than the PN vector.
1019 // This does not cause trouble but is a good indication of a poor
1020 // mesh.
1021 if (skewness > internalSkew)
1022 {
1023 ++nWarnSkew;
1024
1025 if (report)
1026 {
1027 Pout<< "Severe skewness for coupled face " << facei
1028 << " skewness = " << skewness << endl;
1029 }
1030 if (setPtr)
1031 {
1032 setPtr->insert(facei);
1033 }
1034 }
1035
1036 maxSkew = max(maxSkew, skewness);
1037 }
1038 else
1039 {
1041 (
1042 faces,
1043 points,
1044 faceCentres,
1045 faceAreas,
1046
1047 facei,
1048 cellCentres[own[facei]]
1049 );
1050
1051
1052 // Check if the skewness vector is greater than the PN vector.
1053 // This does not cause trouble but is a good indication of a poor
1054 // mesh.
1055 if (skewness > boundarySkew)
1056 {
1057 ++nWarnSkew;
1058
1059 if (report)
1060 {
1061 Pout<< "Severe skewness for boundary face " << facei
1062 << " skewness = " << skewness << endl;
1063 }
1064 if (setPtr)
1065 {
1066 setPtr->insert(facei);
1067 }
1068 }
1069
1070 maxSkew = max(maxSkew, skewness);
1071 }
1072 }
1073
1074 for (const labelPair& baffle : baffles)
1075 {
1076 const label face0 = baffle.first();
1077 const label face1 = baffle.second();
1078
1079 const point& ownCc = cellCentres[own[face0]];
1080 const point& neiCc = cellCentres[own[face1]];
1081
1082 scalar skewness = primitiveMeshTools::faceSkewness
1083 (
1084 faces,
1085 points,
1086 faceCentres,
1087 faceAreas,
1088
1089 face0,
1090 ownCc,
1091 neiCc
1092 );
1093
1094 // Check if the skewness vector is greater than the PN vector.
1095 // This does not cause trouble but is a good indication of a poor
1096 // mesh.
1097 if (skewness > internalSkew)
1098 {
1099 ++nWarnSkew;
1100
1101 if (report)
1102 {
1103 Pout<< "Severe skewness for face " << face0
1104 << " skewness = " << skewness << endl;
1105 }
1106 if (setPtr)
1107 {
1108 setPtr->insert(face0);
1109 }
1110 }
1111
1112 maxSkew = max(maxSkew, skewness);
1113 }
1114
1115
1116 reduce(maxSkew, maxOp<scalar>());
1117 reduce(nWarnSkew, sumOp<label>());
1118
1119 if (nWarnSkew > 0)
1120 {
1121 if (report)
1122 {
1124 << 100*maxSkew
1125 << " percent.\nThis may impair the quality of the result." << nl
1126 << nWarnSkew << " highly skew faces detected."
1127 << endl;
1128 }
1129
1130 return true;
1131 }
1132
1133 if (report)
1134 {
1135 Info<< "Max skewness = " << 100*maxSkew
1136 << " percent. Face skewness OK.\n" << endl;
1137 }
1138
1139 return false;
1140}
1141
1142
1144(
1145 const bool report,
1146 const scalar warnWeight,
1147 const polyMesh& mesh,
1148 const vectorField& cellCentres,
1149 const vectorField& faceCentres,
1150 const vectorField& faceAreas,
1151 const labelList& checkFaces,
1152 const List<labelPair>& baffles,
1153 labelHashSet* setPtr
1154)
1155{
1156 // Warn if the delta factor (0..1) is too large.
1157
1158 const labelList& own = mesh.faceOwner();
1159 const labelList& nei = mesh.faceNeighbour();
1161
1162 // Calculate coupled cell centre
1164
1165 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
1166 {
1167 neiCc[facei-mesh.nInternalFaces()] = cellCentres[own[facei]];
1168 }
1170
1171
1172 scalar minWeight = GREAT;
1173
1174 label nWarnWeight = 0;
1175
1176 for (const label facei : checkFaces)
1177 {
1178 const point& fc = faceCentres[facei];
1179 const vector& fa = faceAreas[facei];
1180
1181 scalar dOwn = mag(fa & (fc-cellCentres[own[facei]]));
1182
1183 if (mesh.isInternalFace(facei))
1184 {
1185 scalar dNei = mag(fa & (cellCentres[nei[facei]]-fc));
1186 scalar weight = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
1187
1188 if (weight < warnWeight)
1189 {
1190 ++nWarnWeight;
1191
1192 if (report)
1193 {
1194 Pout<< "Small weighting factor for face " << facei
1195 << " weight = " << weight << endl;
1196 }
1197 if (setPtr)
1198 {
1199 setPtr->insert(facei);
1200 }
1201 }
1202
1203 minWeight = min(minWeight, weight);
1204 }
1205 else
1206 {
1207 label patchi = patches.whichPatch(facei);
1208
1209 if (patches[patchi].coupled())
1210 {
1211 scalar dNei = mag(fa & (neiCc[facei-mesh.nInternalFaces()]-fc));
1212 scalar weight = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
1213
1214 if (weight < warnWeight)
1215 {
1216 ++nWarnWeight;
1217
1218 if (report)
1219 {
1220 Pout<< "Small weighting factor for face " << facei
1221 << " weight = " << weight << endl;
1222 }
1223 if (setPtr)
1224 {
1225 setPtr->insert(facei);
1226 }
1227 }
1228
1229 minWeight = min(minWeight, weight);
1230 }
1231 }
1232 }
1233
1234 for (const labelPair& baffle : baffles)
1235 {
1236 const label face0 = baffle.first();
1237 const label face1 = baffle.second();
1238
1239 const point& ownCc = cellCentres[own[face0]];
1240 const point& fc = faceCentres[face0];
1241 const vector& fa = faceAreas[face0];
1242
1243 scalar dOwn = mag(fa & (fc-ownCc));
1244 scalar dNei = mag(fa & (cellCentres[own[face1]]-fc));
1245 scalar weight = min(dNei,dOwn)/(dNei+dOwn+VSMALL);
1246
1247 if (weight < warnWeight)
1248 {
1249 ++nWarnWeight;
1250
1251 if (report)
1252 {
1253 Pout<< "Small weighting factor for face " << face0
1254 << " weight = " << weight << endl;
1255 }
1256 if (setPtr)
1257 {
1258 setPtr->insert(face0);
1259 }
1260 }
1261
1262 minWeight = min(minWeight, weight);
1263 }
1264
1265 reduce(minWeight, minOp<scalar>());
1266 reduce(nWarnWeight, sumOp<label>());
1267
1268 if (minWeight < warnWeight)
1269 {
1270 if (report)
1271 {
1273 << minWeight << '.' << nl
1274 << nWarnWeight << " faces with small weights detected."
1275 << endl;
1276 }
1277
1278 return true;
1279 }
1280
1281 if (report)
1282 {
1283 Info<< "Min weight = " << minWeight
1284 << ". Weights OK.\n" << endl;
1285 }
1286
1287 return false;
1288}
1289
1290
1292(
1293 const bool report,
1294 const scalar warnRatio,
1295 const polyMesh& mesh,
1296 const scalarField& cellVolumes,
1297 const labelList& checkFaces,
1298 const List<labelPair>& baffles,
1299 labelHashSet* setPtr
1300)
1301{
1302 // Warn if the volume ratio between neighbouring cells is too large
1303
1304 const labelList& own = mesh.faceOwner();
1305 const labelList& nei = mesh.faceNeighbour();
1307
1308 // Calculate coupled cell vol
1309 scalarField neiVols(mesh.nBoundaryFaces());
1310
1311 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
1312 {
1313 neiVols[facei-mesh.nInternalFaces()] = cellVolumes[own[facei]];
1314 }
1316
1317
1318 scalar minRatio = GREAT;
1319
1320 label nWarnRatio = 0;
1321
1322 for (const label facei : checkFaces)
1323 {
1324 scalar ownVol = mag(cellVolumes[own[facei]]);
1325
1326 scalar neiVol = -GREAT;
1327
1328 if (mesh.isInternalFace(facei))
1329 {
1330 neiVol = mag(cellVolumes[nei[facei]]);
1331 }
1332 else
1333 {
1334 label patchi = patches.whichPatch(facei);
1335
1336 if (patches[patchi].coupled())
1337 {
1338 neiVol = mag(neiVols[facei-mesh.nInternalFaces()]);
1339 }
1340 }
1341
1342 if (neiVol >= 0)
1343 {
1344 scalar ratio = min(ownVol, neiVol) / (max(ownVol, neiVol) + VSMALL);
1345
1346 if (ratio < warnRatio)
1347 {
1348 ++nWarnRatio;
1349
1350 if (report)
1351 {
1352 Pout<< "Small ratio for face " << facei
1353 << " ratio = " << ratio << endl;
1354 }
1355 if (setPtr)
1356 {
1357 setPtr->insert(facei);
1358 }
1359 }
1360
1361 minRatio = min(minRatio, ratio);
1362 }
1363 }
1364
1365 for (const labelPair& baffle : baffles)
1366 {
1367 const label face0 = baffle.first();
1368 const label face1 = baffle.second();
1369
1370 scalar ownVol = mag(cellVolumes[own[face0]]);
1371
1372 scalar neiVol = mag(cellVolumes[own[face1]]);
1373
1374 if (neiVol >= 0)
1375 {
1376 scalar ratio = min(ownVol, neiVol) / (max(ownVol, neiVol) + VSMALL);
1377
1378 if (ratio < warnRatio)
1379 {
1380 ++nWarnRatio;
1381
1382 if (report)
1383 {
1384 Pout<< "Small ratio for face " << face0
1385 << " ratio = " << ratio << endl;
1386 }
1387 if (setPtr)
1388 {
1389 setPtr->insert(face0);
1390 }
1391 }
1392
1393 minRatio = min(minRatio, ratio);
1394 }
1395 }
1396
1397 reduce(minRatio, minOp<scalar>());
1398 reduce(nWarnRatio, sumOp<label>());
1399
1400 if (minRatio < warnRatio)
1401 {
1402 if (report)
1403 {
1405 << minRatio << '.' << nl
1406 << nWarnRatio << " faces with small ratios detected."
1407 << endl;
1408 }
1409
1410 return true;
1411 }
1412
1413 if (report)
1414 {
1415 Info<< "Min ratio = " << minRatio
1416 << ". Ratios OK.\n" << endl;
1417 }
1418
1419 return false;
1420}
1421
1422
1423// Check convexity of angles in a face. Allow a slight non-convexity.
1424// E.g. maxDeg = 10 allows for angles < 190 (or 10 degrees concavity)
1425// (if truly concave and points not visible from face centre the face-pyramid
1426// check in checkMesh will fail)
1428(
1429 const bool report,
1430 const scalar maxDeg,
1431 const polyMesh& mesh,
1432 const vectorField& faceAreas,
1433 const pointField& p,
1434 const labelList& checkFaces,
1435 labelHashSet* setPtr
1436)
1437{
1438 if (maxDeg < -SMALL || maxDeg > 180+SMALL)
1439 {
1441 << "maxDeg should be [0..180] but is now " << maxDeg
1442 << abort(FatalError);
1443 }
1444
1445 const scalar maxSin = Foam::sin(degToRad(maxDeg));
1446
1447 const faceList& fcs = mesh.faces();
1448
1449 scalar maxEdgeSin = 0.0;
1450
1451 label nConcave = 0;
1452
1453 label errorFacei = -1;
1454
1455 for (const label facei : checkFaces)
1456 {
1457 const face& f = fcs[facei];
1458
1459 const vector faceNormal = normalised(faceAreas[facei]);
1460
1461 // Normalized vector from f[size-1] to f[0];
1462 vector ePrev(p[f.first()] - p[f.last()]);
1463 scalar magEPrev = mag(ePrev);
1464 ePrev /= magEPrev + VSMALL;
1465
1466 forAll(f, fp0)
1467 {
1468 // Normalized vector between two consecutive points
1469 vector e10(p[f.nextLabel(fp0)] - p[f.thisLabel(fp0)]);
1470 scalar magE10 = mag(e10);
1471 e10 /= magE10 + VSMALL;
1472
1473 if (magEPrev > SMALL && magE10 > SMALL)
1474 {
1475 vector edgeNormal = ePrev ^ e10;
1476 scalar magEdgeNormal = mag(edgeNormal);
1477
1478 if (magEdgeNormal < maxSin)
1479 {
1480 // Edges (almost) aligned -> face is ok.
1481 }
1482 else
1483 {
1484 // Check normal
1485 edgeNormal /= magEdgeNormal;
1486
1487 if ((edgeNormal & faceNormal) < SMALL)
1488 {
1489 if (facei != errorFacei)
1490 {
1491 // Count only one error per face.
1492 errorFacei = facei;
1493 ++nConcave;
1494 }
1495
1496 maxEdgeSin = max(maxEdgeSin, magEdgeNormal);
1497
1498 if (setPtr)
1499 {
1500 setPtr->insert(facei);
1501 }
1502 }
1503 }
1504 }
1505
1506 ePrev = e10;
1507 magEPrev = magE10;
1508 }
1509 }
1510
1511 reduce(nConcave, sumOp<label>());
1512 reduce(maxEdgeSin, maxOp<scalar>());
1513
1514 if (report)
1515 {
1516 if (maxEdgeSin > SMALL)
1517 {
1518 scalar maxConcaveDegr =
1519 radToDeg(Foam::asin(Foam::min(1.0, maxEdgeSin)));
1520
1521 Info<< "There are " << nConcave
1522 << " faces with concave angles between consecutive"
1523 << " edges. Max concave angle = " << maxConcaveDegr
1524 << " degrees.\n" << endl;
1525 }
1526 else
1527 {
1528 Info<< "All angles in faces are convex or less than " << maxDeg
1529 << " degrees concave.\n" << endl;
1530 }
1531 }
1532
1533 if (nConcave > 0)
1534 {
1535 if (report)
1536 {
1538 << nConcave << " face points with severe concave angle (> "
1539 << maxDeg << " deg) found.\n"
1540 << endl;
1541 }
1542
1543 return true;
1544 }
1545
1546 return false;
1547}
1548
1549
1550// Check twist of faces. Is calculated as the difference between normals of
1551// individual triangles and the cell-cell centre edge.
1553(
1554 const bool report,
1555 const scalar minTwist,
1556 const polyMesh& mesh,
1557 const vectorField& cellCentres,
1558 const vectorField& faceAreas,
1559 const vectorField& faceCentres,
1560 const pointField& p,
1561 const labelList& checkFaces,
1562 labelHashSet* setPtr
1563)
1564{
1565 if (minTwist < -1-SMALL || minTwist > 1+SMALL)
1566 {
1568 << "minTwist should be [-1..1] but is now " << minTwist
1569 << abort(FatalError);
1570 }
1571
1572
1573 const faceList& fcs = mesh.faces();
1574
1575 label nWarped = 0;
1576
1577 const labelList& own = mesh.faceOwner();
1578 const labelList& nei = mesh.faceNeighbour();
1580
1581 // Calculate coupled cell centre
1583
1584 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
1585 {
1586 neiCc[facei-mesh.nInternalFaces()] = cellCentres[own[facei]];
1587 }
1589
1590 for (const label facei : checkFaces)
1591 {
1592 const face& f = fcs[facei];
1593
1594 if (f.size() > 3)
1595 {
1596 vector nf(Zero);
1597
1598 if (mesh.isInternalFace(facei))
1599 {
1600 nf =
1602 (
1603 cellCentres[nei[facei]]
1604 - cellCentres[own[facei]]
1605 );
1606 }
1607 else if (patches[patches.whichPatch(facei)].coupled())
1608 {
1609 nf =
1611 (
1612 neiCc[facei-mesh.nInternalFaces()]
1613 - cellCentres[own[facei]]
1614 );
1615 }
1616 else
1617 {
1618 nf =
1620 (
1621 faceCentres[facei]
1622 - cellCentres[own[facei]]
1623 );
1624 }
1625
1626 if (nf != vector::zero)
1627 {
1628 const point& fc = faceCentres[facei];
1629
1630 forAll(f, fpI)
1631 {
1632 vector triArea
1633 (
1635 (
1636 p[f[fpI]],
1637 p[f.nextLabel(fpI)],
1638 fc
1639 ).areaNormal()
1640 );
1641
1642 scalar magTri = mag(triArea);
1643
1644 if (magTri > VSMALL && ((nf & triArea/magTri) < minTwist))
1645 {
1646 ++nWarped;
1647
1648 if (setPtr)
1649 {
1650 setPtr->insert(facei);
1651 }
1652
1653 break;
1654 }
1655 }
1656 }
1657 }
1658 }
1659
1660 reduce(nWarped, sumOp<label>());
1661
1662 if (report)
1663 {
1664 if (nWarped> 0)
1665 {
1666 Info<< "There are " << nWarped
1667 << " faces with cosine of the angle"
1668 << " between triangle normal and face normal less than "
1669 << minTwist << nl << endl;
1670 }
1671 else
1672 {
1673 Info<< "All faces are flat in that the cosine of the angle"
1674 << " between triangle normal and face normal less than "
1675 << minTwist << nl << endl;
1676 }
1677 }
1678
1679 if (nWarped > 0)
1680 {
1681 if (report)
1682 {
1684 << nWarped << " faces with severe warpage "
1685 << "(cosine of the angle between triangle normal and "
1686 << "face normal < " << minTwist << ") found.\n"
1687 << endl;
1688 }
1689
1690 return true;
1691 }
1692
1693 return false;
1694}
1695
1696
1697// Like checkFaceTwist but compares normals of consecutive triangles.
1699(
1700 const bool report,
1701 const scalar minTwist,
1702 const polyMesh& mesh,
1703 const vectorField& faceAreas,
1704 const vectorField& faceCentres,
1705 const pointField& p,
1706 const labelList& checkFaces,
1707 labelHashSet* setPtr
1708)
1709{
1710 if (minTwist < -1-SMALL || minTwist > 1+SMALL)
1711 {
1713 << "minTwist should be [-1..1] but is now " << minTwist
1714 << abort(FatalError);
1715 }
1716
1717 const faceList& fcs = mesh.faces();
1718
1719 label nWarped = 0;
1720
1721 for (const label facei : checkFaces)
1722 {
1723 const face& f = fcs[facei];
1724
1725 if (f.size() > 3)
1726 {
1727 const point& fc = faceCentres[facei];
1728
1729 // Find starting triangle (at startFp) with non-zero area
1730 label startFp = -1;
1731 vector prevN;
1732
1733 forAll(f, fp)
1734 {
1735 prevN = triPointRef
1736 (
1737 p[f[fp]],
1738 p[f.nextLabel(fp)],
1739 fc
1740 ).areaNormal();
1741
1742 scalar magTri = mag(prevN);
1743
1744 if (magTri > VSMALL)
1745 {
1746 startFp = fp;
1747 prevN /= magTri;
1748 break;
1749 }
1750 }
1751
1752 if (startFp != -1)
1753 {
1754 label fp = startFp;
1755
1756 do
1757 {
1758 fp = f.fcIndex(fp);
1759
1760 vector triN
1761 (
1763 (
1764 p[f[fp]],
1765 p[f.nextLabel(fp)],
1766 fc
1767 ).areaNormal()
1768 );
1769 scalar magTri = mag(triN);
1770
1771 if (magTri > VSMALL)
1772 {
1773 triN /= magTri;
1774
1775 if ((prevN & triN) < minTwist)
1776 {
1777 ++nWarped;
1778
1779 if (setPtr)
1780 {
1781 setPtr->insert(facei);
1782 }
1783
1784 break;
1785 }
1786
1787 prevN = triN;
1788 }
1789 else if (minTwist > 0)
1790 {
1791 ++nWarped;
1792
1793 if (setPtr)
1794 {
1795 setPtr->insert(facei);
1796 }
1797
1798 break;
1799 }
1800 }
1801 while (fp != startFp);
1802 }
1803 }
1804 }
1805
1806
1807 reduce(nWarped, sumOp<label>());
1808
1809 if (report)
1810 {
1811 if (nWarped> 0)
1812 {
1813 Info<< "There are " << nWarped
1814 << " faces with cosine of the angle"
1815 << " between consecutive triangle normals less than "
1816 << minTwist << nl << endl;
1817 }
1818 else
1819 {
1820 Info<< "All faces are flat in that the cosine of the angle"
1821 << " between consecutive triangle normals is less than "
1822 << minTwist << nl << endl;
1823 }
1824 }
1825
1826 if (nWarped > 0)
1827 {
1828 if (report)
1829 {
1831 << nWarped << " faces with severe warpage "
1832 << "(cosine of the angle between consecutive triangle normals"
1833 << " < " << minTwist << ") found.\n"
1834 << endl;
1835 }
1836
1837 return true;
1838 }
1839
1840 return false;
1841}
1842
1843
1845(
1846 const bool report,
1847 const scalar minFlatness,
1848 const polyMesh& mesh,
1849 const vectorField& faceAreas,
1850 const vectorField& faceCentres,
1851 const pointField& p,
1852 const labelList& checkFaces,
1853 labelHashSet* setPtr
1854)
1855{
1856 if (minFlatness < -SMALL || minFlatness > 1+SMALL)
1857 {
1859 << "minFlatness should be [0..1] but is now " << minFlatness
1860 << abort(FatalError);
1861 }
1862
1863 const faceList& fcs = mesh.faces();
1864
1865 label nWarped = 0;
1866
1867 for (const label facei : checkFaces)
1868 {
1869 const face& f = fcs[facei];
1870
1871 if (f.size() > 3)
1872 {
1873 const point& fc = faceCentres[facei];
1874
1875 // Sum triangle areas
1876 scalar sumArea = 0.0;
1877
1878 forAll(f, fp)
1879 {
1880 sumArea += triPointRef
1881 (
1882 p[f[fp]],
1883 p[f.nextLabel(fp)],
1884 fc
1885 ).mag();
1886 }
1887
1888 if (sumArea/mag(faceAreas[facei]) < minFlatness)
1889 {
1890 ++nWarped;
1891 if (setPtr)
1892 {
1893 setPtr->insert(facei);
1894 }
1895 }
1896 }
1897 }
1898
1899 reduce(nWarped, sumOp<label>());
1900
1901 if (report)
1902 {
1903 if (nWarped> 0)
1904 {
1905 Info<< "There are " << nWarped
1906 << " faces with area of individual triangles"
1907 << " compared to overall area less than "
1908 << minFlatness << nl << endl;
1909 }
1910 else
1911 {
1912 Info<< "All faces are flat in that the area of individual triangles"
1913 << " compared to overall area is less than "
1914 << minFlatness << nl << endl;
1915 }
1916 }
1917
1918 if (nWarped > 0)
1919 {
1920 if (report)
1921 {
1923 << nWarped << " non-flat faces "
1924 << "(area of individual triangles"
1925 << " compared to overall area"
1926 << " < " << minFlatness << ") found.\n"
1927 << endl;
1928 }
1929
1930 return true;
1931 }
1932
1933 return false;
1934}
1935
1936
1938(
1939 const bool report,
1940 const scalar minArea,
1941 const polyMesh& mesh,
1942 const vectorField& faceAreas,
1943 const labelList& checkFaces,
1944 labelHashSet* setPtr
1945)
1946{
1947 label nZeroArea = 0;
1948
1949 for (const label facei : checkFaces)
1950 {
1951 if (mag(faceAreas[facei]) < minArea)
1952 {
1953 ++nZeroArea;
1954 if (setPtr)
1955 {
1956 setPtr->insert(facei);
1957 }
1958 }
1959 }
1960
1961
1962 reduce(nZeroArea, sumOp<label>());
1963
1964 if (report)
1965 {
1966 if (nZeroArea > 0)
1967 {
1968 Info<< "There are " << nZeroArea
1969 << " faces with area < " << minArea << '.' << nl << endl;
1970 }
1971 else
1972 {
1973 Info<< "All faces have area > " << minArea << '.' << nl << endl;
1974 }
1975 }
1976
1977 if (nZeroArea > 0)
1978 {
1979 if (report)
1980 {
1982 << nZeroArea << " faces with area < " << minArea
1983 << " found.\n"
1984 << endl;
1985 }
1986
1987 return true;
1988 }
1989
1990 return false;
1991}
1992
1993
1995(
1996 const bool report,
1997 const scalar warnDet,
1998 const polyMesh& mesh,
1999 const vectorField& faceAreas,
2000 const labelList& checkFaces,
2001 const labelList& affectedCells,
2002 labelHashSet* setPtr
2003)
2004{
2005 const cellList& cells = mesh.cells();
2006
2007 scalar minDet = GREAT;
2008 scalar sumDet = 0.0;
2009 label nSumDet = 0;
2010 label nWarnDet = 0;
2011
2012 for (const label celli : affectedCells)
2013 {
2014 const cell& cFaces = cells[celli];
2015
2016 tensor areaSum(Zero);
2017 scalar magAreaSum = 0;
2018
2019 for (const label facei : cFaces)
2020 {
2021 scalar magArea = mag(faceAreas[facei]);
2022
2023 magAreaSum += magArea;
2024 areaSum += faceAreas[facei]*(faceAreas[facei]/(magArea+VSMALL));
2025 }
2026
2027 scalar scaledDet = det(areaSum/(magAreaSum+VSMALL))/0.037037037037037;
2028
2029 minDet = min(minDet, scaledDet);
2030 sumDet += scaledDet;
2031 ++nSumDet;
2032
2033 if (scaledDet < warnDet)
2034 {
2035 ++nWarnDet;
2036 if (setPtr)
2037 {
2038 setPtr->insert(cFaces); // Insert all faces of the cell
2039 }
2040 }
2041 }
2042
2043 reduce(minDet, minOp<scalar>());
2044 reduce(sumDet, sumOp<scalar>());
2045 reduce(nSumDet, sumOp<label>());
2046 reduce(nWarnDet, sumOp<label>());
2047
2048 if (report)
2049 {
2050 if (nSumDet > 0)
2051 {
2052 Info<< "Cell determinant (1 = uniform cube) : average = "
2053 << sumDet / nSumDet << " min = " << minDet << endl;
2054 }
2055
2056 if (nWarnDet > 0)
2057 {
2058 Info<< "There are " << nWarnDet
2059 << " cells with determinant < " << warnDet << '.' << nl
2060 << endl;
2061 }
2062 else
2063 {
2064 Info<< "All faces have determinant > " << warnDet << '.' << nl
2065 << endl;
2066 }
2067 }
2068
2069 if (nWarnDet > 0)
2070 {
2071 if (report)
2072 {
2074 << nWarnDet << " cells with determinant < " << warnDet
2075 << " found.\n"
2076 << endl;
2077 }
2078
2079 return true;
2080 }
2081
2082 return false;
2083}
2084
2085
2087(
2088 const bool report,
2089 const scalar orthWarn,
2090 const labelList& checkFaces,
2091 const List<labelPair>& baffles,
2092 labelHashSet* setPtr
2093) const
2094{
2095 return checkFaceDotProduct
2096 (
2097 report,
2098 orthWarn,
2099 mesh_,
2100 cellCentres_,
2101 faceAreas_,
2102 checkFaces,
2103 baffles,
2104 setPtr
2105 );
2106}
2107
2108
2110(
2111 const bool report,
2112 const scalar minPyrVol,
2113 const pointField& p,
2114 const labelList& checkFaces,
2115 const List<labelPair>& baffles,
2116 labelHashSet* setPtr
2117) const
2118{
2119 return checkFacePyramids
2120 (
2121 report,
2122 minPyrVol,
2123 mesh_,
2124 cellCentres_,
2125 p,
2126 checkFaces,
2127 baffles,
2128 setPtr
2129 );
2130}
2131
2132
2134(
2135 const bool report,
2136 const scalar minTetQuality,
2137 const pointField& p,
2138 const labelList& checkFaces,
2139 const List<labelPair>& baffles,
2140 labelHashSet* setPtr
2141) const
2142{
2143 return checkFaceTets
2144 (
2145 report,
2146 minTetQuality,
2147 mesh_,
2148 cellCentres_,
2149 faceCentres_,
2150 p,
2151 checkFaces,
2152 baffles,
2153 setPtr
2154 );
2155}
2156
2157
2158//bool Foam::polyMeshGeometry::checkFaceSkewness
2159//(
2160// const bool report,
2161// const scalar internalSkew,
2162// const scalar boundarySkew,
2163// const labelList& checkFaces,
2164// const List<labelPair>& baffles,
2165// labelHashSet* setPtr
2166//) const
2167//{
2168// return checkFaceSkewness
2169// (
2170// report,
2171// internalSkew,
2172// boundarySkew,
2173// mesh_,
2174// cellCentres_,
2175// faceCentres_,
2176// faceAreas_,
2177// checkFaces,
2178// baffles,
2179// setPtr
2180// );
2181//}
2182
2183
2185(
2186 const bool report,
2187 const scalar warnWeight,
2188 const labelList& checkFaces,
2189 const List<labelPair>& baffles,
2190 labelHashSet* setPtr
2191) const
2192{
2193 return checkFaceWeights
2194 (
2195 report,
2196 warnWeight,
2197 mesh_,
2198 cellCentres_,
2199 faceCentres_,
2200 faceAreas_,
2201 checkFaces,
2202 baffles,
2203 setPtr
2204 );
2205}
2206
2207
2209(
2210 const bool report,
2211 const scalar warnRatio,
2212 const labelList& checkFaces,
2213 const List<labelPair>& baffles,
2214 labelHashSet* setPtr
2215) const
2216{
2217 return checkVolRatio
2218 (
2219 report,
2220 warnRatio,
2221 mesh_,
2222 cellVolumes_,
2223 checkFaces,
2224 baffles,
2225 setPtr
2226 );
2227}
2228
2229
2231(
2232 const bool report,
2233 const scalar maxDeg,
2234 const pointField& p,
2235 const labelList& checkFaces,
2236 labelHashSet* setPtr
2237) const
2238{
2239 return checkFaceAngles
2240 (
2241 report,
2242 maxDeg,
2243 mesh_,
2244 faceAreas_,
2245 p,
2246 checkFaces,
2247 setPtr
2248 );
2249}
2250
2251
2253(
2254 const bool report,
2255 const scalar minTwist,
2256 const pointField& p,
2257 const labelList& checkFaces,
2258 labelHashSet* setPtr
2259) const
2260{
2261 return checkFaceTwist
2262 (
2263 report,
2264 minTwist,
2265 mesh_,
2266 cellCentres_,
2267 faceAreas_,
2268 faceCentres_,
2269 p,
2270 checkFaces,
2271 setPtr
2272 );
2273}
2274
2275
2277(
2278 const bool report,
2279 const scalar minTwist,
2280 const pointField& p,
2281 const labelList& checkFaces,
2282 labelHashSet* setPtr
2283) const
2284{
2285 return checkTriangleTwist
2286 (
2287 report,
2288 minTwist,
2289 mesh_,
2290 faceAreas_,
2291 faceCentres_,
2292 p,
2293 checkFaces,
2294 setPtr
2295 );
2296}
2297
2298
2300(
2301 const bool report,
2302 const scalar minFlatness,
2303 const pointField& p,
2304 const labelList& checkFaces,
2305 labelHashSet* setPtr
2306) const
2307{
2308 return checkFaceFlatness
2309 (
2310 report,
2311 minFlatness,
2312 mesh_,
2313 faceAreas_,
2314 faceCentres_,
2315 p,
2316 checkFaces,
2317 setPtr
2318 );
2319}
2320
2321
2323(
2324 const bool report,
2325 const scalar minArea,
2326 const labelList& checkFaces,
2327 labelHashSet* setPtr
2328) const
2329{
2330 return checkFaceArea
2331 (
2332 report,
2333 minArea,
2334 mesh_,
2335 faceAreas_,
2336 checkFaces,
2337 setPtr
2338 );
2339}
2340
2341
2343(
2344 const bool report,
2345 const scalar warnDet,
2346 const labelList& checkFaces,
2347 const labelList& affectedCells,
2348 labelHashSet* setPtr
2349) const
2350{
2351 return checkCellDeterminant
2352 (
2353 report,
2354 warnDet,
2355 mesh_,
2356 faceAreas_,
2357 checkFaces,
2358 affectedCells,
2359 setPtr
2360 );
2361}
2362
2363
2364// ************************************************************************* //
label n
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:122
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:216
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
static constexpr direction size() noexcept
The number of elements in the VectorSpace = Ncmpts.
Definition: VectorSpace.H:176
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
tmp< scalarField > faceSkewness() const
Return face skewness.
Definition: cellQuality.C:243
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Updateable mesh geometry and checking routines.
static bool checkFacePyramids(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
See primitiveMesh.
static bool checkTriangleTwist(const bool report, const scalar minTwist, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Consecutive triangle (from face-centre decomposition) normals.
static bool checkFaceTwist(const bool report, const scalar minTwist, const polyMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Triangle (from face-centre decomposition) normal v.s.
static bool checkFaceSkewness(const bool report, const scalar internalSkew, const scalar boundarySkew, const polyMesh &mesh, const pointField &points, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
See primitiveMesh.
void correct()
Take over properties from mesh.
static bool checkFaceArea(const bool report, const scalar minArea, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Small faces.
static bool checkFaceDotProduct(const bool report, const scalar orthWarn, const polyMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
See primitiveMesh.
static bool checkFaceAngles(const bool report, const scalar maxDeg, const polyMesh &mesh, const vectorField &faceAreas, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
See primitiveMesh.
static bool checkFaceWeights(const bool report, const scalar warnWeight, const polyMesh &mesh, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
Interpolation weights (0.5 for regular mesh)
static bool checkVolRatio(const bool report, const scalar warnRatio, const polyMesh &mesh, const scalarField &cellVolumes, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
Cell volume ratio of neighbouring cells (1 for regular mesh)
static bool checkFaceTets(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
See primitiveMesh.
static bool checkFaceFlatness(const bool report, const scalar minFlatness, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Area of faces v.s. sum of triangle areas.
static bool checkCellDeterminant(const bool report, const scalar minDet, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, const labelList &affectedCells, labelHashSet *setPtr)
Area of internal faces v.s. boundary faces.
static labelList affectedCells(const polyMesh &, const labelList &changedFaces)
Helper function: get affected cells from faces.
static label findSharedBasePoint(const polyMesh &mesh, label fI, const point &nCc, scalar tol, bool report=false)
Find the first suitable base point to use for a minimum.
static label findBasePoint(const polyMesh &mesh, label fI, scalar tol, bool report=false)
Find the base point to use for a minimum triangle.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
static scalar boundaryFaceSkewness(const UList< face > &faces, const pointField &p, const vectorField &fCtrs, const vectorField &fAreas, const label facei, const point &ownCc)
Skewness of single boundary face.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
label nInternalFaces() const noexcept
Number of internal faces.
label nFaces() const noexcept
Number of mesh faces.
const cellList & cells() const
static void swapBoundaryFacePositions(const polyMesh &mesh, UList< point > &positions)
Swap coupled positions. Uses eqOp.
Definition: syncTools.H:461
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:445
static void swapBoundaryCellPositions(const polyMesh &mesh, const UList< point > &cellData, List< point > &neighbourCellData)
Swap to obtain neighbour cell positions for all boundary faces.
Definition: syncTools.C:34
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:80
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition: triangleI.H:112
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const polyBoundaryMesh & patches
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
label nPoints
const cellShapeList & cells
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
constexpr scalar pi(M_PI)
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
pyramid< point, const point &, const face & > pyramidPointFaceRef
A pyramid using referred points and faces.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetPointRef.H:47
dimensionedScalar asin(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition: List.H:66
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:71
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Vector< solveScalar > solveVector
Definition: vector.H:64
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59
Unit conversion functions.