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