primitiveMeshGeometry.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 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
30#include "pyramidPointFaceRef.H"
31#include "unitConversion.H"
32#include "triPointRef.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
39}
40
41
42// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43
44void Foam::primitiveMeshGeometry::updateFaceCentresAndAreas
45(
46 const pointField& p,
47 const labelList& changedFaces
48)
49{
50 const faceList& fs = mesh_.faces();
51
52 for (label facei : changedFaces)
53 {
54 const labelList& f = fs[facei];
55 label nPoints = f.size();
56
57 // If the face is a triangle, do a direct calculation for efficiency
58 // and to avoid round-off error-related problems
59 if (nPoints == 3)
60 {
61 faceCentres_[facei] = (1.0/3.0)*(p[f[0]] + p[f[1]] + p[f[2]]);
62 faceAreas_[facei] = 0.5*((p[f[1]] - p[f[0]])^(p[f[2]] - p[f[0]]));
63 }
64 else
65 {
66 vector sumN = Zero;
67 scalar sumA = 0.0;
68 vector sumAc = Zero;
69
70 point fCentre = p[f[0]];
71 for (label pi = 1; pi < nPoints; ++pi)
72 {
73 fCentre += p[f[pi]];
74 }
75
76 fCentre /= nPoints;
77
78 for (label pi = 0; pi < nPoints; ++pi)
79 {
80 const point& nextPoint = p[f[(pi + 1) % nPoints]];
81
82 vector c(p[f[pi]] + nextPoint + fCentre);
83 vector n((nextPoint - p[f[pi]])^(fCentre - p[f[pi]]));
84 scalar a = mag(n);
85
86 sumN += n;
87 sumA += a;
88 sumAc += a*c;
89 }
90
91 faceCentres_[facei] = (1.0/3.0)*sumAc/(sumA + VSMALL);
92 faceAreas_[facei] = 0.5*sumN;
93 }
94 }
95}
96
97
98void Foam::primitiveMeshGeometry::updateCellCentresAndVols
99(
100 const labelList& changedCells,
101 const labelList& changedFaces
102)
103{
104 // Clear the fields for accumulation
105 UIndirectList<vector>(cellCentres_, changedCells) = Zero;
106 UIndirectList<scalar>(cellVolumes_, changedCells) = 0.0;
107
108 const labelList& own = mesh_.faceOwner();
109 const labelList& nei = mesh_.faceNeighbour();
110
111 // first estimate the approximate cell centre as the average of face centres
112
113 vectorField cEst(mesh_.nCells());
114 UIndirectList<vector>(cEst, changedCells) = Zero;
115 scalarField nCellFaces(mesh_.nCells());
116 UIndirectList<scalar>(nCellFaces, changedCells) = 0.0;
117
118 for (label facei : changedFaces)
119 {
120 cEst[own[facei]] += faceCentres_[facei];
121 nCellFaces[own[facei]] += 1;
122
123 if (mesh_.isInternalFace(facei))
124 {
125 cEst[nei[facei]] += faceCentres_[facei];
126 nCellFaces[nei[facei]] += 1;
127 }
128 }
129
130 for (label celli : changedCells)
131 {
132 cEst[celli] /= nCellFaces[celli];
133 }
134
135 for (label facei : changedFaces)
136 {
137 // Calculate 3*face-pyramid volume
138 scalar pyr3Vol = max
139 (
140 faceAreas_[facei] & (faceCentres_[facei] - cEst[own[facei]]),
141 VSMALL
142 );
143
144 // Calculate face-pyramid centre
145 vector pc = (3.0/4.0)*faceCentres_[facei] + (1.0/4.0)*cEst[own[facei]];
146
147 // Accumulate volume-weighted face-pyramid centre
148 cellCentres_[own[facei]] += pyr3Vol*pc;
149
150 // Accumulate face-pyramid volume
151 cellVolumes_[own[facei]] += pyr3Vol;
152
153 if (mesh_.isInternalFace(facei))
154 {
155 // Calculate 3*face-pyramid volume
156 scalar pyr3Vol = max
157 (
158 faceAreas_[facei] & (cEst[nei[facei]] - faceCentres_[facei]),
159 VSMALL
160 );
161
162 // Calculate face-pyramid centre
163 vector pc =
164 (3.0/4.0)*faceCentres_[facei]
165 + (1.0/4.0)*cEst[nei[facei]];
166
167 // Accumulate volume-weighted face-pyramid centre
168 cellCentres_[nei[facei]] += pyr3Vol*pc;
169
170 // Accumulate face-pyramid volume
171 cellVolumes_[nei[facei]] += pyr3Vol;
172 }
173 }
174
175 forAll(changedCells, i)
176 {
177 label celli = changedCells[i];
178
179 cellCentres_[celli] /= cellVolumes_[celli];
180 cellVolumes_[celli] *= (1.0/3.0);
181 }
182}
183
184
186(
187 const labelList& changedFaces
188) const
189{
190 const labelList& own = mesh_.faceOwner();
191 const labelList& nei = mesh_.faceNeighbour();
192
193 labelHashSet affectedCells(2*changedFaces.size());
194
195 for (label facei : changedFaces)
196 {
197 affectedCells.insert(own[facei]);
198
199 if (mesh_.isInternalFace(facei))
200 {
201 affectedCells.insert(nei[facei]);
202 }
203 }
204 return affectedCells.toc();
205}
206
207
208// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
209
211(
212 const primitiveMesh& mesh
213)
214:
215 mesh_(mesh)
216{
217 correct();
218}
219
220
221// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222
224{
225 faceAreas_ = mesh_.faceAreas();
226 faceCentres_ = mesh_.faceCentres();
227 cellCentres_ = mesh_.cellCentres();
228 cellVolumes_ = mesh_.cellVolumes();
229}
230
231
233(
234 const pointField& p,
235 const labelList& changedFaces
236)
237{
238 // Update face quantities
239 updateFaceCentresAndAreas(p, changedFaces);
240 // Update cell quantities from face quantities
241 updateCellCentresAndVols(affectedCells(changedFaces), changedFaces);
242}
243
244
246(
247 const bool report,
248 const scalar orthWarn,
249 const primitiveMesh& mesh,
250 const vectorField& cellCentres,
251 const vectorField& faceAreas,
252 const labelList& checkFaces,
253 labelHashSet* setPtr
254)
255{
256 // for all internal faces check that the d dot S product is positive
257
258 const labelList& own = mesh.faceOwner();
259 const labelList& nei = mesh.faceNeighbour();
260
261 // Severe nonorthogonality threshold
262 const scalar severeNonorthogonalityThreshold = ::cos(degToRad(orthWarn));
263
264 scalar minDDotS = GREAT;
265
266 scalar sumDDotS = 0;
267
268 label severeNonOrth = 0;
269
270 label errorNonOrth = 0;
271
272 forAll(checkFaces, i)
273 {
274 label facei = checkFaces[i];
275
276 if (mesh.isInternalFace(facei))
277 {
278 vector d = cellCentres[nei[facei]] - cellCentres[own[facei]];
279 const vector& s = faceAreas[facei];
280
281 scalar dDotS = (d & s)/(mag(d)*mag(s) + VSMALL);
282
283 if (dDotS < severeNonorthogonalityThreshold)
284 {
285 if (dDotS > SMALL)
286 {
287 if (report)
288 {
289 // Severe non-orthogonality but mesh still OK
290 Pout<< "Severe non-orthogonality for face " << facei
291 << " between cells " << own[facei]
292 << " and " << nei[facei]
293 << ": Angle = " << radToDeg(::acos(dDotS))
294 << " deg." << endl;
295 }
296
297 if (setPtr)
298 {
299 setPtr->insert(facei);
300 }
301
302 ++severeNonOrth;
303 }
304 else
305 {
306 // Non-orthogonality greater than 90 deg
307 if (report)
308 {
310 << "Severe non-orthogonality detected for face "
311 << facei
312 << " between cells " << own[facei] << " and "
313 << nei[facei]
314 << ": Angle = " << radToDeg(::acos(dDotS))
315 << " deg." << endl;
316 }
317
318 ++errorNonOrth;
319
320 if (setPtr)
321 {
322 setPtr->insert(facei);
323 }
324 }
325 }
326
327 if (dDotS < minDDotS)
328 {
329 minDDotS = dDotS;
330 }
331
332 sumDDotS += dDotS;
333 }
334 }
335
336 reduce(minDDotS, minOp<scalar>());
337 reduce(sumDDotS, sumOp<scalar>());
338 reduce(severeNonOrth, sumOp<label>());
339 reduce(errorNonOrth, sumOp<label>());
340
341 label neiSize = nei.size();
342 reduce(neiSize, sumOp<label>());
343
344 // Only report if there are some internal faces
345 if (neiSize > 0)
346 {
347 if (report && minDDotS < severeNonorthogonalityThreshold)
348 {
349 Info<< "Number of non-orthogonality errors: " << errorNonOrth
350 << ". Number of severely non-orthogonal faces: "
351 << severeNonOrth << "." << endl;
352 }
353 }
354
355 if (report)
356 {
357 if (neiSize > 0)
358 {
359 Info<< "Mesh non-orthogonality Max: "
360 << radToDeg(::acos(minDDotS))
361 << " average: " << radToDeg(::acos(sumDDotS/neiSize))
362 << endl;
363 }
364 }
365
366 if (errorNonOrth > 0)
367 {
368 if (report)
369 {
371 << "Error in non-orthogonality detected" << endl;
372 }
373
374 return true;
375 }
376
377 if (report)
378 {
379 Info<< "Non-orthogonality check OK.\n" << endl;
380 }
381
382 return false;
383}
384
385
387(
388 const bool report,
389 const scalar minPyrVol,
390 const primitiveMesh& mesh,
391 const vectorField& cellCentres,
392 const pointField& p,
393 const labelList& checkFaces,
394 labelHashSet* setPtr
395)
396{
397 // check whether face area vector points to the cell with higher label
398 const labelList& own = mesh.faceOwner();
399 const labelList& nei = mesh.faceNeighbour();
400
401 const faceList& f = mesh.faces();
402
403 label nErrorPyrs = 0;
404
405 for (label facei : checkFaces)
406 {
407 // Create the owner pyramid - it will have negative volume
408 scalar pyrVol = pyramidPointFaceRef
409 (
410 f[facei],
411 cellCentres[own[facei]]
412 ).mag(p);
413
414 if (pyrVol > -minPyrVol)
415 {
416 if (report)
417 {
418 Pout<< "bool primitiveMeshGeometry::checkFacePyramids("
419 << "const bool, const scalar, const pointField&"
420 << ", const labelList&, labelHashSet*): "
421 << "face " << facei << " points the wrong way. " << endl
422 << "Pyramid volume: " << -pyrVol
423 << " Face " << f[facei] << " area: " << f[facei].mag(p)
424 << " Owner cell: " << own[facei] << endl
425 << "Owner cell vertex labels: "
426 << mesh.cells()[own[facei]].labels(f)
427 << endl;
428 }
429
430
431 if (setPtr)
432 {
433 setPtr->insert(facei);
434 }
435
436 ++nErrorPyrs;
437 }
438
439 if (mesh.isInternalFace(facei))
440 {
441 // Create the neighbour pyramid - it will have positive volume
442 scalar pyrVol =
443 pyramidPointFaceRef(f[facei], cellCentres[nei[facei]]).mag(p);
444
445 if (pyrVol < minPyrVol)
446 {
447 if (report)
448 {
449 Pout<< "bool primitiveMeshGeometry::checkFacePyramids("
450 << "const bool, const scalar, const pointField&"
451 << ", const labelList&, labelHashSet*): "
452 << "face " << facei << " points the wrong way. " << endl
453 << "Pyramid volume: " << -pyrVol
454 << " Face " << f[facei] << " area: " << f[facei].mag(p)
455 << " Neighbour cell: " << nei[facei] << endl
456 << "Neighbour cell vertex labels: "
457 << mesh.cells()[nei[facei]].labels(f)
458 << endl;
459 }
460
461 if (setPtr)
462 {
463 setPtr->insert(facei);
464 }
465
466 ++nErrorPyrs;
467 }
468 }
469 }
470
471 reduce(nErrorPyrs, sumOp<label>());
472
473 if (nErrorPyrs > 0)
474 {
475 if (report)
476 {
478 << "Error in face pyramids: faces pointing the wrong way!"
479 << endl;
480 }
481
482 return true;
483 }
484
485 if (report)
486 {
487 Info<< "Face pyramids OK.\n" << endl;
488 }
489
490 return false;
491}
492
493
495(
496 const bool report,
497 const scalar internalSkew,
498 const scalar boundarySkew,
499 const primitiveMesh& mesh,
500 const vectorField& cellCentres,
501 const vectorField& faceCentres,
502 const vectorField& faceAreas,
503 const labelList& checkFaces,
504 labelHashSet* setPtr
505)
506{
507 // Warn if the skew correction vector is more than skew times
508 // larger than the face area vector
509
510 const labelList& own = mesh.faceOwner();
511 const labelList& nei = mesh.faceNeighbour();
512
513 scalar maxSkew = 0;
514
515 label nWarnSkew = 0;
516
517 for (label facei : checkFaces)
518 {
519 if (mesh.isInternalFace(facei))
520 {
521 scalar dOwn = mag(faceCentres[facei] - cellCentres[own[facei]]);
522 scalar dNei = mag(faceCentres[facei] - cellCentres[nei[facei]]);
523
524 point faceIntersection =
525 cellCentres[own[facei]]*dNei/(dOwn+dNei)
526 + cellCentres[nei[facei]]*dOwn/(dOwn+dNei);
527
528 scalar skewness =
529 mag(faceCentres[facei] - faceIntersection)
530 / (
531 mag(cellCentres[nei[facei]]-cellCentres[own[facei]])
532 + VSMALL
533 );
534
535 // Check if the skewness vector is greater than the PN vector.
536 // This does not cause trouble but is a good indication of a poor
537 // mesh.
538 if (skewness > internalSkew)
539 {
540 if (report)
541 {
542 Pout<< "Severe skewness for face " << facei
543 << " skewness = " << skewness << endl;
544 }
545
546 if (setPtr)
547 {
548 setPtr->insert(facei);
549 }
550
551 ++nWarnSkew;
552 }
553
554 if (skewness > maxSkew)
555 {
556 maxSkew = skewness;
557 }
558 }
559 else
560 {
561 // Boundary faces: consider them to have only skewness error.
562 // (i.e. treat as if mirror cell on other side)
563
564 const vector faceNormal = normalised(faceAreas[facei]);
565
566 vector dOwn = faceCentres[facei] - cellCentres[own[facei]];
567
568 vector dWall = faceNormal*(faceNormal & dOwn);
569
570 point faceIntersection = cellCentres[own[facei]] + dWall;
571
572 scalar skewness =
573 mag(faceCentres[facei] - faceIntersection)
574 /(2*mag(dWall) + VSMALL);
575
576 // Check if the skewness vector is greater than the PN vector.
577 // This does not cause trouble but is a good indication of a poor
578 // mesh.
579 if (skewness > boundarySkew)
580 {
581 if (report)
582 {
583 Pout<< "Severe skewness for boundary face " << facei
584 << " skewness = " << skewness << endl;
585 }
586
587 if (setPtr)
588 {
589 setPtr->insert(facei);
590 }
591
592 ++nWarnSkew;
593 }
594
595 if (skewness > maxSkew)
596 {
597 maxSkew = skewness;
598 }
599 }
600 }
601
602 reduce(maxSkew, maxOp<scalar>());
603 reduce(nWarnSkew, sumOp<label>());
604
605 if (nWarnSkew > 0)
606 {
607 if (report)
608 {
610 << 100*maxSkew
611 << " percent.\nThis may impair the quality of the result." << nl
612 << nWarnSkew << " highly skew faces detected."
613 << endl;
614 }
615
616 return true;
617 }
618
619 if (report)
620 {
621 Info<< "Max skewness = " << 100*maxSkew
622 << " percent. Face skewness OK.\n" << endl;
623 }
624
625 return false;
626}
627
628
630(
631 const bool report,
632 const scalar warnWeight,
633 const primitiveMesh& mesh,
634 const vectorField& cellCentres,
635 const vectorField& faceCentres,
636 const vectorField& faceAreas,
637 const labelList& checkFaces,
638 labelHashSet* setPtr
639)
640{
641 // Warn if the delta factor (0..1) is too large.
642
643 const labelList& own = mesh.faceOwner();
644 const labelList& nei = mesh.faceNeighbour();
645
646 scalar minWeight = GREAT;
647
648 label nWarnWeight = 0;
649
650 for (label facei : checkFaces)
651 {
652 if (mesh.isInternalFace(facei))
653 {
654 const point& fc = faceCentres[facei];
655
656 scalar dOwn = mag(faceAreas[facei] & (fc-cellCentres[own[facei]]));
657 scalar dNei = mag(faceAreas[facei] & (cellCentres[nei[facei]]-fc));
658
659 scalar weight = min(dNei,dOwn)/(dNei+dOwn);
660
661 if (weight < warnWeight)
662 {
663 if (report)
664 {
665 Pout<< "Small weighting factor for face " << facei
666 << " weight = " << weight << endl;
667 }
668
669 if (setPtr)
670 {
671 setPtr->insert(facei);
672 }
673
674 ++nWarnWeight;
675 }
676
677 minWeight = min(minWeight, weight);
678 }
679 }
680
681 reduce(minWeight, minOp<scalar>());
682 reduce(nWarnWeight, sumOp<label>());
683
684 if (minWeight < warnWeight)
685 {
686 if (report)
687 {
689 << minWeight << '.' << nl
690 << nWarnWeight << " faces with small weights detected."
691 << endl;
692 }
693
694 return true;
695 }
696
697 if (report)
698 {
699 Info<< "Min weight = " << minWeight
700 << " percent. Weights OK.\n" << endl;
701 }
702
703 return false;
704}
705
706
707// Check convexity of angles in a face. Allow a slight non-convexity.
708// E.g. maxDeg = 10 allows for angles < 190 (or 10 degrees concavity)
709// (if truly concave and points not visible from face centre the face-pyramid
710// check in checkMesh will fail)
712(
713 const bool report,
714 const scalar maxDeg,
715 const primitiveMesh& mesh,
716 const vectorField& faceAreas,
717 const pointField& p,
718 const labelList& checkFaces,
719 labelHashSet* setPtr
720)
721{
722 if (maxDeg < -SMALL || maxDeg > 180+SMALL)
723 {
725 << "maxDeg should be [0..180] but is now " << maxDeg
726 << abort(FatalError);
727 }
728
729 const scalar maxSin = Foam::sin(degToRad(maxDeg));
730
731 const faceList& fcs = mesh.faces();
732
733 scalar maxEdgeSin = 0.0;
734
735 label nConcave = 0;
736
737 label errorFacei = -1;
738
739 for (label facei : checkFaces)
740 {
741 const face& f = fcs[facei];
742
743 const vector faceNormal = normalised(faceAreas[facei]);
744
745 // Normalized vector from f[size-1] to f[0];
746 vector ePrev(p[f.first()] - p[f.last()]);
747 scalar magEPrev = mag(ePrev);
748 ePrev /= magEPrev + VSMALL;
749
750 forAll(f, fp0)
751 {
752 // Normalized vector between two consecutive points
753 vector e10(p[f.nextLabel(fp0)] - p[f.thisLabel(fp0)]);
754 scalar magE10 = mag(e10);
755 e10 /= magE10 + VSMALL;
756
757 if (magEPrev > SMALL && magE10 > SMALL)
758 {
759 vector edgeNormal = ePrev ^ e10;
760 scalar magEdgeNormal = mag(edgeNormal);
761
762 if (magEdgeNormal < maxSin)
763 {
764 // Edges (almost) aligned -> face is ok.
765 }
766 else
767 {
768 // Check normal
769 edgeNormal /= magEdgeNormal;
770
771 if ((edgeNormal & faceNormal) < SMALL)
772 {
773 if (facei != errorFacei)
774 {
775 // Count only one error per face.
776 errorFacei = facei;
777 ++nConcave;
778 }
779
780 if (setPtr)
781 {
782 setPtr->insert(facei);
783 }
784
785 maxEdgeSin = max(maxEdgeSin, magEdgeNormal);
786 }
787 }
788 }
789
790 ePrev = e10;
791 magEPrev = magE10;
792 }
793 }
794
795 reduce(nConcave, sumOp<label>());
796 reduce(maxEdgeSin, maxOp<scalar>());
797
798 if (report)
799 {
800 if (maxEdgeSin > SMALL)
801 {
802 scalar maxConcaveDegr =
803 radToDeg(Foam::asin(Foam::min(1.0, maxEdgeSin)));
804
805 Info<< "There are " << nConcave
806 << " faces with concave angles between consecutive"
807 << " edges. Max concave angle = " << maxConcaveDegr
808 << " degrees.\n" << endl;
809 }
810 else
811 {
812 Info<< "All angles in faces are convex or less than " << maxDeg
813 << " degrees concave.\n" << endl;
814 }
815 }
816
817 if (nConcave > 0)
818 {
819 if (report)
820 {
822 << nConcave << " face points with severe concave angle (> "
823 << maxDeg << " deg) found.\n"
824 << endl;
825 }
826
827 return true;
828 }
829
830 return false;
831}
832
833
837//bool Foam::primitiveMeshGeometry::checkFaceFlatness
838//(
839// const bool report,
840// const scalar warnFlatness,
841// const primitiveMesh& mesh,
842// const vectorField& faceAreas,
843// const vectorField& faceCentres,
844// const pointField& p,
845// const labelList& checkFaces,
846// labelHashSet* setPtr
847//)
848//{
849// if (warnFlatness < 0 || warnFlatness > 1)
850// {
851// FatalErrorInFunction
852// << "warnFlatness should be [0..1] but is now " << warnFlatness
853// << abort(FatalError);
854// }
855//
856//
857// const faceList& fcs = mesh.faces();
858//
859// // Areas are calculated as the sum of areas. (see
860// // primitiveMeshFaceCentresAndAreas.C)
861//
862// label nWarped = 0;
863//
864// scalar minFlatness = GREAT;
865// scalar sumFlatness = 0;
866// label nSummed = 0;
867//
868// forAll(checkFaces, i)
869// {
870// label facei = checkFaces[i];
871//
872// const face& f = fcs[facei];
873//
874// scalar magArea = mag(faceAreas[facei]);
875//
876// if (f.size() > 3 && magArea > VSMALL)
877// {
878// const point& fc = faceCentres[facei];
879//
880// // Calculate the sum of magnitude of areas and compare to
881// // magnitude of sum of areas.
882//
883// scalar sumA = 0.0;
884//
885// forAll(f, fp)
886// {
887// const point& thisPoint = p[f[fp]];
888// const point& nextPoint = p[f.nextLabel(fp)];
889//
890// // Triangle around fc.
891// vector n = 0.5*((nextPoint - thisPoint)^(fc - thisPoint));
892// sumA += mag(n);
893// }
894//
895// scalar flatness = magArea / (sumA+VSMALL);
896//
897// sumFlatness += flatness;
898// ++nSummed;
899//
900// minFlatness = min(minFlatness, flatness);
901//
902// if (flatness < warnFlatness)
903// {
904// ++nWarped;
905//
906// if (setPtr)
907// {
908// setPtr->insert(facei);
909// }
910// }
911// }
912// }
913//
914//
915// reduce(nWarped, sumOp<label>());
916// reduce(minFlatness, minOp<scalar>());
917//
918// reduce(nSummed, sumOp<label>());
919// reduce(sumFlatness, sumOp<scalar>());
920//
921// if (report)
922// {
923// if (nSummed > 0)
924// {
925// Info<< "Face flatness (1 = flat, 0 = butterfly) : average = "
926// << sumFlatness / nSummed << " min = " << minFlatness << endl;
927// }
928//
929// if (nWarped> 0)
930// {
931// Info<< "There are " << nWarped
932// << " faces with ratio between projected and actual area < "
933// << warnFlatness
934// << ".\nMinimum ratio (minimum flatness, maximum warpage) = "
935// << minFlatness << nl << endl;
936// }
937// else
938// {
939// Info<< "All faces are flat in that the ratio between projected"
940// << " and actual area is > " << warnFlatness << nl << endl;
941// }
942// }
943//
944// if (nWarped > 0)
945// {
946// if (report)
947// {
948// WarningInFunction
949// << nWarped << " faces with severe warpage (flatness < "
950// << warnFlatness << ") found.\n"
951// << endl;
952// }
953//
954// return true;
955// }
956//
957// return false;
958//}
959
960
961// Check twist of faces. Is calculated as the difference between areas of
962// individual triangles and the overall area of the face (which ifself is
963// is the average of the areas of the individual triangles).
965(
966 const bool report,
967 const scalar minTwist,
968 const primitiveMesh& mesh,
969 const vectorField& faceAreas,
970 const vectorField& faceCentres,
971 const pointField& p,
972 const labelList& checkFaces,
973 labelHashSet* setPtr
974)
975{
976 if (minTwist < -1-SMALL || minTwist > 1+SMALL)
977 {
979 << "minTwist should be [-1..1] but is now " << minTwist
980 << abort(FatalError);
981 }
982
983
984 const faceList& fcs = mesh.faces();
985
986 // Areas are calculated as the sum of areas. (see
987 // primitiveMeshFaceCentresAndAreas.C)
988
989 label nWarped = 0;
990
991 for (const label facei : checkFaces)
992 {
993 const face& f = fcs[facei];
994
995 const scalar magArea = mag(faceAreas[facei]);
996
997 if (f.size() > 3 && magArea > VSMALL)
998 {
999 const vector nf = faceAreas[facei] / magArea;
1000
1001 const point& fc = faceCentres[facei];
1002
1003 forAll(f, fpI)
1004 {
1005 const vector triArea
1006 (
1008 (
1009 p[f[fpI]],
1010 p[f.nextLabel(fpI)],
1011 fc
1012 ).areaNormal()
1013 );
1014
1015 const scalar magTri = mag(triArea);
1016
1017 if (magTri > VSMALL && ((nf & triArea/magTri) < minTwist))
1018 {
1019 ++nWarped;
1020
1021 if (setPtr)
1022 {
1023 setPtr->insert(facei);
1024 }
1025 }
1026 }
1027 }
1028 }
1029
1030
1031 reduce(nWarped, sumOp<label>());
1032
1033 if (report)
1034 {
1035 if (nWarped> 0)
1036 {
1037 Info<< "There are " << nWarped
1038 << " faces with cosine of the angle"
1039 << " between triangle normal and face normal less than "
1040 << minTwist << nl << endl;
1041 }
1042 else
1043 {
1044 Info<< "All faces are flat in that the cosine of the angle"
1045 << " between triangle normal and face normal less than "
1046 << minTwist << nl << endl;
1047 }
1048 }
1049
1050 if (nWarped > 0)
1051 {
1052 if (report)
1053 {
1055 << nWarped << " faces with severe warpage "
1056 << "(cosine of the angle between triangle normal and "
1057 << "face normal < " << minTwist << ") found.\n"
1058 << endl;
1059 }
1060
1061 return true;
1062 }
1063
1064 return false;
1065}
1066
1067
1069(
1070 const bool report,
1071 const scalar minArea,
1072 const primitiveMesh& mesh,
1073 const vectorField& faceAreas,
1074 const labelList& checkFaces,
1075 labelHashSet* setPtr
1076)
1077{
1078 label nZeroArea = 0;
1079
1080 for (label facei : checkFaces)
1081 {
1082 if (mag(faceAreas[facei]) < minArea)
1083 {
1084 if (setPtr)
1085 {
1086 setPtr->insert(facei);
1087 }
1088 ++nZeroArea;
1089 }
1090 }
1091
1092
1093 reduce(nZeroArea, sumOp<label>());
1094
1095 if (report)
1096 {
1097 if (nZeroArea > 0)
1098 {
1099 Info<< "There are " << nZeroArea
1100 << " faces with area < " << minArea << '.' << nl << endl;
1101 }
1102 else
1103 {
1104 Info<< "All faces have area > " << minArea << '.' << nl << endl;
1105 }
1106 }
1107
1108 if (nZeroArea > 0)
1109 {
1110 if (report)
1111 {
1113 << nZeroArea << " faces with area < " << minArea
1114 << " found.\n"
1115 << endl;
1116 }
1117
1118 return true;
1119 }
1120
1121 return false;
1122}
1123
1124
1126(
1127 const bool report,
1128 const scalar warnDet,
1129 const primitiveMesh& mesh,
1130 const vectorField& faceAreas,
1131 const labelList& checkFaces,
1132 const labelList& affectedCells,
1133 labelHashSet* setPtr
1134)
1135{
1136 const cellList& cells = mesh.cells();
1137
1138 scalar minDet = GREAT;
1139 scalar sumDet = 0.0;
1140 label nSumDet = 0;
1141 label nWarnDet = 0;
1142
1143 for (label celli : affectedCells)
1144 {
1145 const cell& cFaces = cells[celli];
1146
1147 tensor areaSum(Zero);
1148 scalar magAreaSum = 0;
1149
1150 forAll(cFaces, cFacei)
1151 {
1152 label facei = cFaces[cFacei];
1153
1154 scalar magArea = mag(faceAreas[facei]);
1155
1156 magAreaSum += magArea;
1157 areaSum += faceAreas[facei]*(faceAreas[facei]/magArea);
1158 }
1159
1160 scalar scaledDet = det(areaSum/magAreaSum)/0.037037037037037;
1161
1162 minDet = min(minDet, scaledDet);
1163 sumDet += scaledDet;
1164 ++nSumDet;
1165
1166 if (scaledDet < warnDet)
1167 {
1168 if (setPtr)
1169 {
1170 // Insert all faces of the cell.
1171 forAll(cFaces, cFacei)
1172 {
1173 label facei = cFaces[cFacei];
1174 setPtr->insert(facei);
1175 }
1176 }
1177 ++nWarnDet;
1178 }
1179 }
1180
1181 reduce(minDet, minOp<scalar>());
1182 reduce(sumDet, sumOp<scalar>());
1183 reduce(nSumDet, sumOp<label>());
1184 reduce(nWarnDet, sumOp<label>());
1185
1186 if (report)
1187 {
1188 if (nSumDet > 0)
1189 {
1190 Info<< "Cell determinant (1 = uniform cube) : average = "
1191 << sumDet / nSumDet << " min = " << minDet << endl;
1192 }
1193
1194 if (nWarnDet > 0)
1195 {
1196 Info<< "There are " << nWarnDet
1197 << " cells with determinant < " << warnDet << '.' << nl
1198 << endl;
1199 }
1200 else
1201 {
1202 Info<< "All faces have determinant > " << warnDet << '.' << nl
1203 << endl;
1204 }
1205 }
1206
1207 if (nWarnDet > 0)
1208 {
1209 if (report)
1210 {
1212 << nWarnDet << " cells with determinant < " << warnDet
1213 << " found.\n"
1214 << endl;
1215 }
1216
1217 return true;
1218 }
1219
1220 return false;
1221}
1222
1223
1225(
1226 const bool report,
1227 const scalar orthWarn,
1228 const labelList& checkFaces,
1229 labelHashSet* setPtr
1230) const
1231{
1232 return checkFaceDotProduct
1233 (
1234 report,
1235 orthWarn,
1236 mesh_,
1237 cellCentres_,
1238 faceAreas_,
1239 checkFaces,
1240 setPtr
1241 );
1242}
1243
1244
1246(
1247 const bool report,
1248 const scalar minPyrVol,
1249 const pointField& p,
1250 const labelList& checkFaces,
1251 labelHashSet* setPtr
1252) const
1253{
1254 return checkFacePyramids
1255 (
1256 report,
1257 minPyrVol,
1258 mesh_,
1259 cellCentres_,
1260 p,
1261 checkFaces,
1262 setPtr
1263 );
1264}
1265
1266
1268(
1269 const bool report,
1270 const scalar internalSkew,
1271 const scalar boundarySkew,
1272 const labelList& checkFaces,
1273 labelHashSet* setPtr
1274) const
1275{
1276 return checkFaceSkewness
1277 (
1278 report,
1279 internalSkew,
1280 boundarySkew,
1281 mesh_,
1282 cellCentres_,
1283 faceCentres_,
1284 faceAreas_,
1285 checkFaces,
1286 setPtr
1287 );
1288}
1289
1290
1292(
1293 const bool report,
1294 const scalar warnWeight,
1295 const labelList& checkFaces,
1296 labelHashSet* setPtr
1297) const
1298{
1299 return checkFaceWeights
1300 (
1301 report,
1302 warnWeight,
1303 mesh_,
1304 cellCentres_,
1305 faceCentres_,
1306 faceAreas_,
1307 checkFaces,
1308 setPtr
1309 );
1310}
1311
1312
1314(
1315 const bool report,
1316 const scalar maxDeg,
1317 const pointField& p,
1318 const labelList& checkFaces,
1319 labelHashSet* setPtr
1320) const
1321{
1322 return checkFaceAngles
1323 (
1324 report,
1325 maxDeg,
1326 mesh_,
1327 faceAreas_,
1328 p,
1329 checkFaces,
1330 setPtr
1331 );
1332}
1333
1334
1335//bool Foam::primitiveMeshGeometry::checkFaceFlatness
1336//(
1337// const bool report,
1338// const scalar warnFlatness,
1339// const pointField& p,
1340// const labelList& checkFaces,
1341// labelHashSet* setPtr
1342//) const
1343//{
1344// return checkFaceFlatness
1345// (
1346// report,
1347// warnFlatness,
1348// mesh_,
1349// faceAreas_,
1350// faceCentres_,
1351// p,
1352// checkFaces,
1353// setPtr
1354// );
1355//}
1356
1357
1359(
1360 const bool report,
1361 const scalar minTwist,
1362 const pointField& p,
1363 const labelList& checkFaces,
1364 labelHashSet* setPtr
1365) const
1366{
1367 return checkFaceTwist
1368 (
1369 report,
1370 minTwist,
1371 mesh_,
1372 faceAreas_,
1373 faceCentres_,
1374 p,
1375 checkFaces,
1376 setPtr
1377 );
1378}
1379
1380
1382(
1383 const bool report,
1384 const scalar minArea,
1385 const labelList& checkFaces,
1386 labelHashSet* setPtr
1387) const
1388{
1389 return checkFaceArea
1390 (
1391 report,
1392 minArea,
1393 mesh_,
1394 faceAreas_,
1395 checkFaces,
1396 setPtr
1397 );
1398}
1399
1400
1402(
1403 const bool report,
1404 const scalar warnDet,
1405 const labelList& checkFaces,
1406 const labelList& affectedCells,
1407 labelHashSet* setPtr
1408) const
1409{
1410 return checkCellDeterminant
1411 (
1412 report,
1413 warnDet,
1414 mesh_,
1415 faceAreas_,
1416 checkFaces,
1417 affectedCells,
1418 setPtr
1419 );
1420}
1421
1422
1423// ************************************************************************* //
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
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
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
Updateable mesh geometry + checking routines.
static bool checkFaceTwist(const bool report, const scalar minTwist, const primitiveMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
static bool checkFaceAngles(const bool report, const scalar maxDeg, const primitiveMesh &mesh, const vectorField &faceAreas, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
static bool checkFaceSkewness(const bool report, const scalar internalSkew, const scalar boundarySkew, const primitiveMesh &mesh, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
void correct()
Take over properties from mesh.
static bool checkFaceWeights(const bool report, const scalar warnWeight, const primitiveMesh &mesh, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
labelList affectedCells(const labelList &changedFaces) const
Helper function: get affected cells from faces.
static bool checkCellDeterminant(const bool report, const scalar minDet, const primitiveMesh &, const vectorField &faceAreas, const labelList &checkFaces, const labelList &affectedCells, labelHashSet *setPtr)
static bool checkFaceArea(const bool report, const scalar minArea, const primitiveMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
static bool checkFaceDotProduct(const bool report, const scalar orthWarn, const primitiveMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
static bool checkFacePyramids(const bool report, const scalar minPyrVol, const primitiveMesh &, const vectorField &cellCentres, const pointField &p, const labelList &checkFaces, labelHashSet *)
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
virtual const faceList & faces() const =0
Return faces.
const cellList & cells() const
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
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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
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
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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)
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
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
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
Unit conversion functions.