Loading...
Searching...
No Matches
polyMeshCheck.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) 2012-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2020 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 "polyMesh.H"
30#include "polyMeshTools.H"
31#include "unitConversion.H"
32#include "syncTools.H"
33
34// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35
36bool Foam::polyMesh::checkFaceOrthogonality
37(
38 const vectorField& fAreas,
39 const vectorField& cellCtrs,
40 const bool report,
41 const bool detailedReport,
42 labelHashSet* setPtr
43) const
44{
45 DebugInFunction << "Checking mesh non-orthogonality" << endl;
46
47 const labelList& own = faceOwner();
48 const labelList& nei = faceNeighbour();
49
50 // Calculate orthogonality for all internal and coupled boundary faces
51 // (1 for uncoupled boundary faces)
52 tmp<scalarField> tortho = polyMeshTools::faceOrthogonality
53 (
54 *this,
55 fAreas,
56 cellCtrs
57 );
58 const scalarField& ortho = tortho.ref();
59
60 // Severe nonorthogonality threshold
61 const scalar severeNonorthogonalityThreshold =
63
64
65 scalar minDDotS = GREAT;
66 scalar sumDDotS = 0.0;
67 label nSummed = 0;
68 label severeNonOrth = 0;
69 label errorNonOrth = 0;
70
71
72 // Statistics only for internal and masters of coupled faces
73 bitSet isMasterFace(syncTools::getInternalOrMasterFaces(*this));
74
75 forAll(ortho, facei)
76 {
77 if (ortho[facei] < severeNonorthogonalityThreshold)
78 {
79 if (ortho[facei] > SMALL)
80 {
81 if (setPtr)
82 {
83 setPtr->insert(facei);
84 }
85
86 severeNonOrth++;
87 }
88 else
89 {
90 // Error : non-ortho too large
91 if (setPtr)
92 {
93 setPtr->insert(facei);
94 }
95 if (detailedReport && errorNonOrth == 0)
96 {
97 // Non-orthogonality greater than 90 deg
99 << "Severe non-orthogonality for face "
100 << facei
101 << " between cells " << own[facei]
102 << " and " << nei[facei]
103 << ": Angle = "
104 << radToDeg(::acos(clamp(ortho[facei], -1, 1)))
105 << " deg." << endl;
106 }
107
108 errorNonOrth++;
109 }
110 }
111
112 if (isMasterFace.test(facei))
113 {
114 minDDotS = min(minDDotS, ortho[facei]);
115 sumDDotS += ortho[facei];
116 nSummed++;
117 }
118 }
119
120 reduce(minDDotS, minOp<scalar>());
121 reduce(sumDDotS, sumOp<scalar>());
122 reduce(nSummed, sumOp<label>());
123 reduce(severeNonOrth, sumOp<label>());
124 reduce(errorNonOrth, sumOp<label>());
125
126 dictionary& meshDict = const_cast<dictionary&>(data().meshDict());
127
128 if (nSummed > 0)
129 {
130 scalar maxNonOrth = radToDeg(::acos(clamp(minDDotS, -1, 1)));
131 scalar aveNonOrth = radToDeg(::acos(clamp(sumDDotS/nSummed, -1, 1)));
132
133 meshDict.set("maxNonOrth", maxNonOrth);
134 meshDict.set("aveNonOrth", aveNonOrth);
135
136 if (debug || report)
137 {
138 Info<< " Mesh non-orthogonality Max: " << maxNonOrth
139 << " average: " << aveNonOrth
140 << endl;
141 }
142 }
143
144 if (severeNonOrth > 0)
145 {
146 meshDict.set("thresholdNonOrth", primitiveMesh::nonOrthThreshold_);
147 meshDict.set("severeNonOrth", severeNonOrth);
148
149 if (debug || report)
150 {
151 Info<< " *Number of severely non-orthogonal (> "
152 << primitiveMesh::nonOrthThreshold_ << " degrees) faces: "
153 << severeNonOrth << "." << endl;
154 }
155 }
156
157 if (errorNonOrth > 0)
158 {
159 meshDict.set("errorNonOrth", errorNonOrth);
160
161 if (debug || report)
162 {
163 Info<< " ***Number of non-orthogonality errors: "
164 << errorNonOrth << "." << endl;
165 }
166
167 return true;
168 }
169
170 if (debug || report)
171 {
172 Info<< " Non-orthogonality check OK." << endl;
173 }
174
175 return false;
176}
177
178
179bool Foam::polyMesh::checkFaceSkewness
180(
181 const pointField& points,
182 const vectorField& fCtrs,
183 const vectorField& fAreas,
184 const vectorField& cellCtrs,
185 const bool report,
186 const bool detailedReport,
187 labelHashSet* setPtr
188) const
189{
190 DebugInFunction << "Checking face skewness" << endl;
191
192 const labelList& own = faceOwner();
193 const labelList& nei = faceNeighbour();
194
195 // Warn if the skew correction vector is more than skewWarning times
196 // larger than the face area vector
197
199 (
200 *this,
201 points,
202 fCtrs,
203 fAreas,
204 cellCtrs
205 );
206 const scalarField& skew = tskew.ref();
207
208 scalar maxSkew = max(skew);
209 label nWarnSkew = 0;
210
211 // Statistics only for all faces except slave coupled faces
212 bitSet isMasterFace(syncTools::getMasterFaces(*this));
213
214 forAll(skew, facei)
215 {
216 // Check if the skewness vector is greater than the PN vector.
217 // This does not cause trouble but is a good indication of a poor mesh.
218 if (skew[facei] > skewThreshold_)
219 {
220 if (setPtr)
221 {
222 setPtr->insert(facei);
223 }
224 if (detailedReport && nWarnSkew == 0)
225 {
226 // Non-orthogonality greater than 90 deg
227 if (isInternalFace(facei))
228 {
230 << "Severe skewness " << skew[facei]
231 << " for face " << facei
232 << " between cells " << own[facei]
233 << " and " << nei[facei];
234 }
235 else
236 {
238 << "Severe skewness " << skew[facei]
239 << " for boundary face " << facei
240 << " on cell " << own[facei];
241 }
242 }
243
244 if (isMasterFace.test(facei))
245 {
246 ++nWarnSkew;
247 }
248 }
249 }
250
251 reduce(maxSkew, maxOp<scalar>());
252 reduce(nWarnSkew, sumOp<label>());
253
254 dictionary& meshDict = const_cast<dictionary&>(data().meshDict());
255 meshDict.set("maxSkew", maxSkew);
256
257 if (nWarnSkew > 0)
258 {
259 meshDict.set("nWarnSkew", nWarnSkew);
260
261 if (debug || report)
262 {
263 Info<< " ***Max skewness = " << maxSkew
264 << ", " << nWarnSkew << " highly skew faces detected"
265 " which may impair the quality of the results"
266 << endl;
267 }
268
269 return true;
270 }
271
272 if (debug || report)
273 {
274 Info<< " Max skewness = " << maxSkew << " OK." << endl;
275 }
276
277 return false;
278}
279
280
281bool Foam::polyMesh::checkEdgeAlignment
282(
283 const pointField& p,
284 const bool report,
286 labelHashSet* setPtr
287) const
288{
289 // Check 1D/2Dness of edges. Gets passed the non-empty directions and
290 // checks all edges in the mesh whether they:
291 // - have no component in a non-empty direction or
292 // - are only in a single non-empty direction.
293 // Empty direction info is passed in as a vector of labels (synchronised)
294 // which are 1 if the direction is non-empty, 0 if it is.
295
296 DebugInFunction << "Checking edge alignment" << endl;
297
298 label nDirs = 0;
299 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
300 {
301 if (directions[cmpt] == 1)
302 {
303 nDirs++;
304 }
305 else if (directions[cmpt] != 0)
306 {
308 << "directions should contain 0 or 1 but is now " << directions
309 << exit(FatalError);
310 }
311 }
312
313 if (nDirs == vector::nComponents)
314 {
315 return false;
316 }
317
318
319 const faceList& fcs = faces();
320
321 EdgeMap<label> edgesInError;
322
323 forAll(fcs, facei)
324 {
325 const face& f = fcs[facei];
326
327 forAll(f, fp)
328 {
329 label p0 = f[fp];
330 label p1 = f.nextLabel(fp);
331 if (p0 < p1)
332 {
333 vector d(p[p1]-p[p0]);
334 scalar magD = mag(d);
335
336 if (magD > ROOTVSMALL)
337 {
338 d /= magD;
339
340 // Check how many empty directions are used by the edge.
341 label nEmptyDirs = 0;
342 label nNonEmptyDirs = 0;
343 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
344 {
345 if (mag(d[cmpt]) > 1e-6)
346 {
347 if (directions[cmpt] == 0)
348 {
349 nEmptyDirs++;
350 }
351 else
352 {
353 nNonEmptyDirs++;
354 }
355 }
356 }
357
358 if (nEmptyDirs == 0)
359 {
360 // Purely in ok directions.
361 }
362 else if (nEmptyDirs == 1)
363 {
364 // Ok if purely in empty directions.
365 if (nNonEmptyDirs > 0)
366 {
367 edgesInError.insert(edge(p0, p1), facei);
368 }
369 }
370 else if (nEmptyDirs > 1)
371 {
372 // Always an error
373 edgesInError.insert(edge(p0, p1), facei);
374 }
375 }
376 }
377 }
378 }
379
380 label nErrorEdges = returnReduce(edgesInError.size(), sumOp<label>());
381
382 if (nErrorEdges > 0)
383 {
384 dictionary& meshDict = const_cast<dictionary&>(data().meshDict());
385 meshDict.set("nErrorAlignedEdges", nErrorEdges);
386
387 if (debug || report)
388 {
389 Info<< " ***Number of edges not aligned with or perpendicular to "
390 << "non-empty directions: " << nErrorEdges << endl;
391 }
392
393 if (setPtr)
394 {
395 setPtr->reserve(2*edgesInError.size());
396 forAllConstIters(edgesInError, iter)
397 {
398 setPtr->insert(iter.key().first());
399 setPtr->insert(iter.key().second());
400 }
401 }
402
403 return true;
404 }
405
406 if (debug || report)
407 {
408 Info<< " All edges aligned with or perpendicular to "
409 << "non-empty directions." << endl;
410 }
411
412 return false;
413}
414
415
416bool Foam::polyMesh::checkCellDeterminant
417(
418 const vectorField& faceAreas,
419 const bool report,
420 labelHashSet* setPtr,
421 const Vector<label>& meshD
422) const
423{
424 const scalar warnDet = 1e-3;
425
426 DebugInFunction << "Checking for under-determined cells" << endl;
427
429 (
430 *this,
431 meshD,
432 faceAreas,
434 );
435 scalarField& cellDeterminant = tcellDeterminant.ref();
436
437
438 label nErrorCells = 0;
439 scalar minDet = min(cellDeterminant);
440 scalar sumDet = sum(cellDeterminant);
441
442 forAll(cellDeterminant, celli)
443 {
444 if (cellDeterminant[celli] < warnDet)
445 {
446 if (setPtr)
447 {
448 setPtr->insert(celli);
449 }
450
451 nErrorCells++;
452 }
453 }
454
455 reduce(nErrorCells, sumOp<label>());
456 reduce(minDet, minOp<scalar>());
457 reduce(sumDet, sumOp<scalar>());
458 label nSummed = returnReduce(cellDeterminant.size(), sumOp<label>());
459
460 dictionary& meshDict = const_cast<dictionary&>(data().meshDict());
461
462 if (nSummed > 0)
463 {
464 meshDict.set("minDeterminant", minDet);
465 meshDict.set("aveDeterminant", sumDet/nSummed);
466
467 if (debug || report)
468 {
469 Info<< " Cell determinant (wellposedness) : minimum: " << minDet
470 << " average: " << sumDet/nSummed
471 << endl;
472 }
473 }
474
475 if (nErrorCells > 0)
476 {
477 meshDict.set("thresholdDeterminant", warnDet);
478 meshDict.set("nErrorDeterminant", nErrorCells);
479
480 if (debug || report)
481 {
482 Info<< " ***Cells with small determinant (< "
483 << warnDet << ") found, number of cells: "
484 << nErrorCells << endl;
485 }
486
487 return true;
488 }
489
490 if (debug || report)
491 {
492 Info<< " Cell determinant check OK." << endl;
493 }
494
495 return false;
496}
497
498
499bool Foam::polyMesh::checkFaceWeight
500(
501 const vectorField& fCtrs,
502 const vectorField& fAreas,
503 const vectorField& cellCtrs,
504 const bool report,
505 const scalar minWeight,
506 labelHashSet* setPtr
507) const
508{
509 DebugInFunction << "Checking for low face interpolation weights" << endl;
510
512 (
513 *this,
514 fCtrs,
515 fAreas,
516 cellCtrs
517 );
518 scalarField& faceWght = tfaceWght.ref();
519
520
521 label nErrorFaces = 0;
522 scalar minDet = GREAT;
523 scalar sumDet = 0.0;
524 label nSummed = 0;
525
526 // Statistics only for internal and masters of coupled faces
527 bitSet isMasterFace(syncTools::getInternalOrMasterFaces(*this));
528
529 forAll(faceWght, facei)
530 {
531 if (faceWght[facei] < minWeight)
532 {
533 // Note: insert both sides of coupled faces
534 if (setPtr)
535 {
536 setPtr->insert(facei);
537 }
538
539 nErrorFaces++;
540 }
541
542 // Note: statistics only on master of coupled faces
543 if (isMasterFace.test(facei))
544 {
545 minDet = min(minDet, faceWght[facei]);
546 sumDet += faceWght[facei];
547 nSummed++;
548 }
549 }
550
551 reduce(nErrorFaces, sumOp<label>());
552 reduce(minDet, minOp<scalar>());
553 reduce(sumDet, sumOp<scalar>());
554 reduce(nSummed, sumOp<label>());
555
556 dictionary& meshDict = const_cast<dictionary&>(data().meshDict());
557
558 if (nSummed > 0)
559 {
560 meshDict.set("minFaceWeight", minDet);
561 meshDict.set("aveFaceWeight", sumDet/nSummed);
562
563 if (debug || report)
564 {
565 Info<< " Face interpolation weight : minimum: " << minDet
566 << " average: " << sumDet/nSummed
567 << endl;
568 }
569 }
570
571 if (nErrorFaces > 0)
572 {
573 meshDict.set("thresholdFaceWeight", minWeight);
574 meshDict.set("nErrorFaceWeight", nErrorFaces);
575
576 if (debug || report)
577 {
578 Info<< " ***Faces with small interpolation weight (< " << minWeight
579 << ") found, number of faces: "
580 << nErrorFaces << endl;
581 }
582
583 return true;
584 }
585
586 if (debug || report)
587 {
588 Info<< " Face interpolation weight check OK." << endl;
589 }
590
591 return false;
592}
593
594
595bool Foam::polyMesh::checkVolRatio
596(
597 const scalarField& cellVols,
598 const bool report,
599 const scalar minRatio,
600 labelHashSet* setPtr
601) const
602{
603 DebugInFunction << "Checking for volume ratio < " << minRatio << endl;
604
606 scalarField& volRatio = tvolRatio.ref();
607
608
609 label nErrorFaces = 0;
610 scalar minDet = GREAT;
611 scalar sumDet = 0.0;
612 label nSummed = 0;
613
614 // Statistics only for internal and masters of coupled faces
615 bitSet isMasterFace(syncTools::getInternalOrMasterFaces(*this));
616
617 forAll(volRatio, facei)
618 {
619 if (volRatio[facei] < minRatio)
620 {
621 // Note: insert both sides of coupled faces
622 if (setPtr)
623 {
624 setPtr->insert(facei);
625 }
626
627 nErrorFaces++;
628 }
629
630 // Note: statistics only on master of coupled faces
631 if (isMasterFace.test(facei))
632 {
633 minDet = min(minDet, volRatio[facei]);
634 sumDet += volRatio[facei];
635 nSummed++;
636 }
637 }
638
639 reduce(nErrorFaces, sumOp<label>());
640 reduce(minDet, minOp<scalar>());
641 reduce(sumDet, sumOp<scalar>());
642 reduce(nSummed, sumOp<label>());
643
644 dictionary& meshDict = const_cast<dictionary&>(data().meshDict());
645
646 if (nSummed > 0)
647 {
648 meshDict.set("minFaceVolumeRatio", minDet);
649 meshDict.set("aveFaceVolumeRatio", sumDet/nSummed);
650
651 if (debug || report)
652 {
653 Info<< " Face volume ratio : minimum: " << minDet
654 << " average: " << sumDet/nSummed
655 << endl;
656 }
657 }
658
659 if (nErrorFaces > 0)
660 {
661 meshDict.set("thresholdFaceVolumeRatio", minRatio);
662 meshDict.set("nErrorFaceVolumeRatio", nErrorFaces);
663
664 if (debug || report)
665 {
666 Info<< " ***Faces with small volume ratio (< " << minRatio
667 << ") found, number of faces: "
668 << nErrorFaces << endl;
669 }
670
671 return true;
672 }
673
674 if (debug || report)
675 {
676 Info<< " Face volume ratio check OK." << endl;
677 }
678
679 return false;
680}
681
682
683bool Foam::polyMesh::checkFaceOrthogonality
684(
685 const bool report,
686 labelHashSet* setPtr
687) const
688{
689 return checkFaceOrthogonality
690 (
691 faceAreas(),
692 cellCentres(),
693 report,
694 false, // detailedReport
695 setPtr
696 );
697}
698
699
700bool Foam::polyMesh::checkFaceSkewness
701(
702 const bool report,
703 labelHashSet* setPtr
704) const
705{
706 return checkFaceSkewness
707 (
708 points(),
709 faceCentres(),
710 faceAreas(),
711 cellCentres(),
712 report,
713 false, // detailedReport
714 setPtr
715 );
716}
717
718
719bool Foam::polyMesh::checkEdgeAlignment
720(
721 const bool report,
723 labelHashSet* setPtr
724) const
725{
726 return checkEdgeAlignment
727 (
728 points(),
729 report,
731 setPtr
732 );
733}
734
735
736bool Foam::polyMesh::checkCellDeterminant
737(
738 const bool report,
739 labelHashSet* setPtr
740) const
741{
742 return checkCellDeterminant
743 (
744 faceAreas(),
745 report,
746 setPtr,
747 geometricD()
748 );
749}
750
751
752bool Foam::polyMesh::checkFaceWeight
753(
754 const bool report,
755 const scalar minWeight,
756 labelHashSet* setPtr
757) const
758{
759 return checkFaceWeight
760 (
761 faceCentres(),
762 faceAreas(),
763 cellCentres(),
764 report,
765 minWeight,
766 setPtr
767 );
768}
769
770
771bool Foam::polyMesh::checkVolRatio
772(
773 const bool report,
774 const scalar minRatio,
776) const
777{
778 return checkVolRatio(cellVolumes(), report, minRatio, setPtr);
779}
780
781
783(
784 const pointField& newPoints,
785 const bool report,
786 const bool detailedReport
787) const
788{
789 if (debug || report)
790 {
791 Pout<< "bool polyMesh::checkMeshMotion("
792 << "const pointField&, const bool, const bool) const: "
793 << "checking mesh motion" << endl;
794 }
795
796 vectorField fCtrs(nFaces());
797 vectorField fAreas(nFaces());
798
800 (
801 *this,
802 newPoints,
803 fCtrs,
804 fAreas
805 );
806
807 // Check cell volumes and calculate new cell centres
808 vectorField cellCtrs(nCells());
809 scalarField cellVols(nCells());
810
812 (
813 *this,
814 fCtrs,
815 fAreas,
816 cellCtrs,
818 );
819
820 // Check cell volumes
821 bool error = checkCellVolumes
822 (
823 cellVols, // vols
824 report, // report
825 detailedReport, // detailedReport
826 nullptr // setPtr
827 );
828
829
830 // Check face areas
831 bool areaError = checkFaceAreas
832 (
833 fAreas,
834 report, // report
835 detailedReport, // detailedReport,
836 nullptr // setPtr
837 );
838 error = error || areaError;
839
840
841 // Check pyramid volumes
842 bool pyrVolError = checkFacePyramids
843 (
844 newPoints,
845 cellCtrs,
846 report, // report,
847 detailedReport, // detailedReport,
848 -SMALL, // minPyrVol
849 nullptr // setPtr
850 );
851 error = error || pyrVolError;
852
853
854 // Check face non-orthogonality
855 bool nonOrthoError = checkFaceOrthogonality
856 (
857 fAreas,
858 cellCtrs,
859 report, // report
860 detailedReport, // detailedReport
861 nullptr // setPtr
862 );
863 error = error || nonOrthoError;
864
865
866 if (!error && (debug || report))
867 {
868 Pout<< "Mesh motion check OK." << endl;
869 }
870
871 return error;
872}
873
874
875// ************************************************************************* //
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
static constexpr direction nComponents
Number of components in this vector space.
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition directions.H:67
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition error.H:74
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
dictionary & meshDict()
Return the dictionary of mesh data, typically populated by the polyMesh::checkXXX functions,...
Definition meshState.C:113
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate orthogonality field. (1 for fully orthogonal, < 1 for non-orthogonal).
static tmp< scalarField > faceSkewness(const polyMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
static tmp< scalarField > faceWeights(const polyMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate interpolation factors field.
static tmp< scalarField > volRatio(const polyMesh &mesh, const scalarField &vol)
Generate volume ratio field.
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
virtual const meshState & data() const noexcept
Const reference to the mesh and solver state data.
Definition polyMesh.H:559
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
virtual bool checkMeshMotion(const pointField &newPoints, const bool report=false, const bool detailedReport=false) const
Check mesh motion for correctness given motion points.
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition polyMesh.C:850
static void makeFaceCentresAndAreas(const UList< face > &faces, const pointField &p, vectorField &fCtrs, vectorField &fAreas)
Calculate face centres and areas for specified faces.
static tmp< scalarField > cellDeterminant(const primitiveMesh &mesh, const Vector< label > &directions, const vectorField &faceAreas, const bitSet &internalOrCoupledFace)
Generate cell determinant field. Normalised to 1 for an internal cube.
static void makeCellCentresAndVols(const primitiveMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols)
Calculate cell centres and volumes from face properties.
const vectorField & faceCentres() const
const scalarField & cellVolumes() const
const vectorField & cellCentres() const
bool checkFacePyramids(const pointField &points, const vectorField &ctrs, const bool report, const bool detailedReport, const scalar minPyrVol, labelHashSet *setPtr) const
Check face pyramid volume.
static scalar nonOrthThreshold_
Non-orthogonality warning threshold in deg.
bool checkCellVolumes(const scalarField &vols, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative cell volumes.
bool checkFaceAreas(const vectorField &faceAreas, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative face areas.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const vectorField & faceAreas() const
static bitSet getMasterFaces(const polyMesh &mesh)
Get per face whether it is uncoupled or a master of a coupled set of faces.
Definition syncTools.C:119
static bitSet getInternalOrCoupledFaces(const polyMesh &mesh)
Get per face whether it is internal or coupled.
Definition syncTools.C:165
static bitSet getInternalOrMasterFaces(const polyMesh &mesh)
Get per face whether it is internal or a master of a coupled set of faces.
Definition syncTools.C:139
A class for managing temporary objects.
Definition tmp.H:75
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const IOdictionary & meshDict
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
List< face > faceList
List of faces.
Definition faceListFwd.H:41
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.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedTensor skew(const dimensionedTensor &dt)
dimensionedScalar acos(const dimensionedScalar &ds)
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
const scalarField & cellVols
Unit conversion functions.