polyMeshCheck.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "polyMesh.H"
27 #include "polyMeshTools.H"
28 #include "unitConversion.H"
29 #include "syncTools.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 bool Foam::polyMesh::checkFaceOrthogonality
34 (
35  const vectorField& fAreas,
36  const vectorField& cellCtrs,
37  const bool report,
38  const bool detailedReport,
39  labelHashSet* setPtr
40 ) const
41 {
42  if (debug)
43  {
44  InfoInFunction << "Checking mesh non-orthogonality" << endl;
45  }
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  PackedBoolList 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(min(1.0, max(-1.0, ortho[facei]))))
105  << " deg." << endl;
106  }
107 
108  errorNonOrth++;
109  }
110  }
111 
112  if (isMasterFace[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  if (debug || report)
127  {
128  if (nSummed > 0)
129  {
130  if (debug || report)
131  {
132  Info<< " Mesh non-orthogonality Max: "
133  << radToDeg(::acos(min(1.0, max(-1.0, minDDotS))))
134  << " average: "
135  << radToDeg(::acos(min(1.0, max(-1.0, sumDDotS/nSummed))))
136  << endl;
137  }
138  }
139 
140  if (severeNonOrth > 0)
141  {
142  Info<< " *Number of severely non-orthogonal (> "
143  << primitiveMesh::nonOrthThreshold_ << " degrees) faces: "
144  << severeNonOrth << "." << endl;
145  }
146  }
147 
148  if (errorNonOrth > 0)
149  {
150  if (debug || report)
151  {
152  Info<< " ***Number of non-orthogonality errors: "
153  << errorNonOrth << "." << endl;
154  }
155 
156  return true;
157  }
158  else
159  {
160  if (debug || report)
161  {
162  Info<< " Non-orthogonality check OK." << endl;
163  }
164 
165  return false;
166  }
167 }
168 
169 
170 bool Foam::polyMesh::checkFaceSkewness
171 (
172  const pointField& points,
173  const vectorField& fCtrs,
174  const vectorField& fAreas,
175  const vectorField& cellCtrs,
176  const bool report,
177  const bool detailedReport,
178  labelHashSet* setPtr
179 ) const
180 {
181  if (debug)
182  {
183  InfoInFunction << "Checking face skewness" << endl;
184  }
185 
186  const labelList& own = faceOwner();
187  const labelList& nei = faceNeighbour();
188 
189  // Warn if the skew correction vector is more than skewWarning times
190  // larger than the face area vector
191 
192  tmp<scalarField> tskew = polyMeshTools::faceSkewness
193  (
194  *this,
195  points,
196  fCtrs,
197  fAreas,
198  cellCtrs
199  );
200  const scalarField& skew = tskew.ref();
201 
202  scalar maxSkew = max(skew);
203  label nWarnSkew = 0;
204 
205  // Statistics only for all faces except slave coupled faces
206  PackedBoolList isMasterFace(syncTools::getMasterFaces(*this));
207 
208  forAll(skew, facei)
209  {
210  // Check if the skewness vector is greater than the PN vector.
211  // This does not cause trouble but is a good indication of a poor mesh.
212  if (skew[facei] > skewThreshold_)
213  {
214  if (setPtr)
215  {
216  setPtr->insert(facei);
217  }
218  if (detailedReport && nWarnSkew == 0)
219  {
220  // Non-orthogonality greater than 90 deg
221  if (isInternalFace(facei))
222  {
224  << "Severe skewness " << skew[facei]
225  << " for face " << facei
226  << " between cells " << own[facei]
227  << " and " << nei[facei];
228  }
229  else
230  {
232  << "Severe skewness " << skew[facei]
233  << " for boundary face " << facei
234  << " on cell " << own[facei];
235  }
236  }
237 
238  if (isMasterFace[facei])
239  {
240  nWarnSkew++;
241  }
242  }
243  }
244 
245  reduce(maxSkew, maxOp<scalar>());
246  reduce(nWarnSkew, sumOp<label>());
247 
248  if (nWarnSkew > 0)
249  {
250  if (debug || report)
251  {
252  Info<< " ***Max skewness = " << maxSkew
253  << ", " << nWarnSkew << " highly skew faces detected"
254  " which may impair the quality of the results"
255  << endl;
256  }
257 
258  return true;
259  }
260  else
261  {
262  if (debug || report)
263  {
264  Info<< " Max skewness = " << maxSkew << " OK." << endl;
265  }
266 
267  return false;
268  }
269 }
270 
271 
272 bool Foam::polyMesh::checkEdgeAlignment
273 (
274  const pointField& p,
275  const bool report,
276  const Vector<label>& directions,
277  labelHashSet* setPtr
278 ) const
279 {
280  // Check 1D/2Dness of edges. Gets passed the non-empty directions and
281  // checks all edges in the mesh whether they:
282  // - have no component in a non-empty direction or
283  // - are only in a single non-empty direction.
284  // Empty direction info is passed in as a vector of labels (synchronised)
285  // which are 1 if the direction is non-empty, 0 if it is.
286 
287  if (debug)
288  {
289  InfoInFunction << "Checking edge alignment" << endl;
290  }
291 
292  label nDirs = 0;
293  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
294  {
295  if (directions[cmpt] == 1)
296  {
297  nDirs++;
298  }
299  else if (directions[cmpt] != 0)
300  {
302  << "directions should contain 0 or 1 but is now " << directions
303  << exit(FatalError);
304  }
305  }
306 
307  if (nDirs == vector::nComponents)
308  {
309  return false;
310  }
311 
312 
313  const faceList& fcs = faces();
314 
315  EdgeMap<label> edgesInError;
316 
317  forAll(fcs, facei)
318  {
319  const face& f = fcs[facei];
320 
321  forAll(f, fp)
322  {
323  label p0 = f[fp];
324  label p1 = f.nextLabel(fp);
325  if (p0 < p1)
326  {
327  vector d(p[p1]-p[p0]);
328  scalar magD = mag(d);
329 
330  if (magD > rootVSmall)
331  {
332  d /= magD;
333 
334  // Check how many empty directions are used by the edge.
335  label nEmptyDirs = 0;
336  label nNonEmptyDirs = 0;
337  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
338  {
339  if (mag(d[cmpt]) > 1e-6)
340  {
341  if (directions[cmpt] == 0)
342  {
343  nEmptyDirs++;
344  }
345  else
346  {
347  nNonEmptyDirs++;
348  }
349  }
350  }
351 
352  if (nEmptyDirs == 0)
353  {
354  // Purely in ok directions.
355  }
356  else if (nEmptyDirs == 1)
357  {
358  // Ok if purely in empty directions.
359  if (nNonEmptyDirs > 0)
360  {
361  edgesInError.insert(edge(p0, p1), facei);
362  }
363  }
364  else if (nEmptyDirs > 1)
365  {
366  // Always an error
367  edgesInError.insert(edge(p0, p1), facei);
368  }
369  }
370  }
371  }
372  }
373 
374  label nErrorEdges = returnReduce(edgesInError.size(), sumOp<label>());
375 
376  if (nErrorEdges > 0)
377  {
378  if (debug || report)
379  {
380  Info<< " ***Number of edges not aligned with or perpendicular to "
381  << "non-empty directions: " << nErrorEdges << endl;
382  }
383 
384  if (setPtr)
385  {
386  setPtr->resize(2*edgesInError.size());
387  forAllConstIter(EdgeMap<label>, edgesInError, iter)
388  {
389  setPtr->insert(iter.key()[0]);
390  setPtr->insert(iter.key()[1]);
391  }
392  }
393 
394  return true;
395  }
396  else
397  {
398  if (debug || report)
399  {
400  Info<< " All edges aligned with or perpendicular to "
401  << "non-empty directions." << endl;
402  }
403  return false;
404  }
405 }
406 
407 
408 bool Foam::polyMesh::checkCellDeterminant
409 (
410  const vectorField& faceAreas,
411  const bool report,
412  labelHashSet* setPtr,
413  const Vector<label>& meshD
414 ) const
415 {
416  const scalar warnDet = 1e-3;
417 
418  if (debug)
419  {
420  InfoInFunction << "Checking for under-determined cells" << endl;
421  }
422 
423  tmp<scalarField> tcellDeterminant = primitiveMeshTools::cellDeterminant
424  (
425  *this,
426  meshD,
427  faceAreas,
429  );
430  scalarField& cellDeterminant = tcellDeterminant.ref();
431 
432 
433  label nErrorCells = 0;
434  scalar minDet = min(cellDeterminant);
435  scalar sumDet = sum(cellDeterminant);
436 
437  forAll(cellDeterminant, celli)
438  {
439  if (cellDeterminant[celli] < warnDet)
440  {
441  if (setPtr)
442  {
443  setPtr->insert(celli);
444  }
445 
446  nErrorCells++;
447  }
448  }
449 
450  reduce(nErrorCells, sumOp<label>());
451  reduce(minDet, minOp<scalar>());
452  reduce(sumDet, sumOp<scalar>());
453  label nSummed = returnReduce(cellDeterminant.size(), sumOp<label>());
454 
455  if (debug || report)
456  {
457  if (nSummed > 0)
458  {
459  Info<< " Cell determinant (wellposedness) : minimum: " << minDet
460  << " average: " << sumDet/nSummed
461  << endl;
462  }
463  }
464 
465  if (nErrorCells > 0)
466  {
467  if (debug || report)
468  {
469  Info<< " ***Cells with small determinant (< "
470  << warnDet << ") found, number of cells: "
471  << nErrorCells << endl;
472  }
473 
474  return true;
475  }
476  else
477  {
478  if (debug || report)
479  {
480  Info<< " Cell determinant check OK." << endl;
481  }
482 
483  return false;
484  }
485 
486  return false;
487 }
488 
489 
490 bool Foam::polyMesh::checkFaceWeight
491 (
492  const vectorField& fCtrs,
493  const vectorField& fAreas,
494  const vectorField& cellCtrs,
495  const bool report,
496  const scalar minWeight,
497  labelHashSet* setPtr
498 ) const
499 {
500  if (debug)
501  {
502  InfoInFunction << "Checking for low face interpolation weights" << endl;
503  }
504 
505  tmp<scalarField> tfaceWght = polyMeshTools::faceWeights
506  (
507  *this,
508  fCtrs,
509  fAreas,
510  cellCtrs
511  );
512  scalarField& faceWght = tfaceWght.ref();
513 
514 
515  label nErrorFaces = 0;
516  scalar minDet = great;
517  scalar sumDet = 0.0;
518  label nSummed = 0;
519 
520  // Statistics only for internal and masters of coupled faces
521  PackedBoolList isMasterFace(syncTools::getInternalOrMasterFaces(*this));
522 
523  forAll(faceWght, facei)
524  {
525  if (faceWght[facei] < minWeight)
526  {
527  // Note: insert both sides of coupled faces
528  if (setPtr)
529  {
530  setPtr->insert(facei);
531  }
532 
533  nErrorFaces++;
534  }
535 
536  // Note: statistics only on master of coupled faces
537  if (isMasterFace[facei])
538  {
539  minDet = min(minDet, faceWght[facei]);
540  sumDet += faceWght[facei];
541  nSummed++;
542  }
543  }
544 
545  reduce(nErrorFaces, sumOp<label>());
546  reduce(minDet, minOp<scalar>());
547  reduce(sumDet, sumOp<scalar>());
548  reduce(nSummed, sumOp<label>());
549 
550  if (debug || report)
551  {
552  if (nSummed > 0)
553  {
554  Info<< " Face interpolation weight : minimum: " << minDet
555  << " average: " << sumDet/nSummed
556  << endl;
557  }
558  }
559 
560  if (nErrorFaces > 0)
561  {
562  if (debug || report)
563  {
564  Info<< " ***Faces with small interpolation weight (< " << minWeight
565  << ") found, number of faces: "
566  << nErrorFaces << endl;
567  }
568 
569  return true;
570  }
571  else
572  {
573  if (debug || report)
574  {
575  Info<< " Face interpolation weight check OK." << endl;
576  }
577 
578  return false;
579  }
580 
581  return false;
582 }
583 
584 
585 bool Foam::polyMesh::checkVolRatio
586 (
587  const scalarField& cellVols,
588  const bool report,
589  const scalar minRatio,
590  labelHashSet* setPtr
591 ) const
592 {
593  if (debug)
594  {
595  InfoInFunction << "Checking for volume ratio < " << minRatio << endl;
596  }
597 
598  tmp<scalarField> tvolRatio = polyMeshTools::volRatio(*this, cellVols);
599  scalarField& volRatio = tvolRatio.ref();
600 
601 
602  label nErrorFaces = 0;
603  scalar minDet = great;
604  scalar sumDet = 0.0;
605  label nSummed = 0;
606 
607  // Statistics only for internal and masters of coupled faces
608  PackedBoolList isMasterFace(syncTools::getInternalOrMasterFaces(*this));
609 
610  forAll(volRatio, facei)
611  {
612  if (volRatio[facei] < minRatio)
613  {
614  // Note: insert both sides of coupled faces
615  if (setPtr)
616  {
617  setPtr->insert(facei);
618  }
619 
620  nErrorFaces++;
621  }
622 
623  // Note: statistics only on master of coupled faces
624  if (isMasterFace[facei])
625  {
626  minDet = min(minDet, volRatio[facei]);
627  sumDet += volRatio[facei];
628  nSummed++;
629  }
630  }
631 
632  reduce(nErrorFaces, sumOp<label>());
633  reduce(minDet, minOp<scalar>());
634  reduce(sumDet, sumOp<scalar>());
635  reduce(nSummed, sumOp<label>());
636 
637  if (debug || report)
638  {
639  if (nSummed > 0)
640  {
641  Info<< " Face volume ratio : minimum: " << minDet
642  << " average: " << sumDet/nSummed
643  << endl;
644  }
645  }
646 
647  if (nErrorFaces > 0)
648  {
649  if (debug || report)
650  {
651  Info<< " ***Faces with small volume ratio (< " << minRatio
652  << ") found, number of faces: "
653  << nErrorFaces << endl;
654  }
655 
656  return true;
657  }
658  else
659  {
660  if (debug || report)
661  {
662  Info<< " Face volume ratio check OK." << endl;
663  }
664 
665  return false;
666  }
667 
668  return false;
669 }
670 
671 
672 bool Foam::polyMesh::checkFaceOrthogonality
673 (
674  const bool report,
675  labelHashSet* setPtr
676 ) const
677 {
678  return checkFaceOrthogonality
679  (
680  faceAreas(),
681  cellCentres(),
682  report,
683  false, // detailedReport
684  setPtr
685  );
686 }
687 
688 
689 bool Foam::polyMesh::checkFaceSkewness
690 (
691  const bool report,
692  labelHashSet* setPtr
693 ) const
694 {
695  return checkFaceSkewness
696  (
697  points(),
698  faceCentres(),
699  faceAreas(),
700  cellCentres(),
701  report,
702  false, // detailedReport
703  setPtr
704  );
705 }
706 
707 
708 bool Foam::polyMesh::checkEdgeAlignment
709 (
710  const bool report,
711  const Vector<label>& directions,
712  labelHashSet* setPtr
713 ) const
714 {
715  return checkEdgeAlignment
716  (
717  points(),
718  report,
719  directions,
720  setPtr
721  );
722 }
723 
724 
725 bool Foam::polyMesh::checkCellDeterminant
726 (
727  const bool report,
728  labelHashSet* setPtr
729 ) const
730 {
731  return checkCellDeterminant
732  (
733  faceAreas(),
734  report,
735  setPtr,
736  geometricD()
737  );
738 }
739 
740 
741 bool Foam::polyMesh::checkFaceWeight
742 (
743  const bool report,
744  const scalar minWeight,
745  labelHashSet* setPtr
746 ) const
747 {
748  return checkFaceWeight
749  (
750  faceCentres(),
751  faceAreas(),
752  cellCentres(),
753  report,
754  minWeight,
755  setPtr
756  );
757 }
758 
759 
760 bool Foam::polyMesh::checkVolRatio
761 (
762  const bool report,
763  const scalar minRatio,
764  labelHashSet* setPtr
765 ) const
766 {
767  return checkVolRatio(cellVolumes(), report, minRatio, setPtr);
768 }
769 
770 
772 (
773  const pointField& newPoints,
774  const bool report,
775  const bool detailedReport
776 ) const
777 {
778  if (debug || report)
779  {
780  Pout<< "bool polyMesh::checkMeshMotion("
781  << "const pointField&, const bool, const bool) const: "
782  << "checking mesh motion" << endl;
783  }
784 
785  vectorField fCtrs(nFaces());
786  vectorField fAreas(nFaces());
787  scalarField magfAreas(nFaces());
788 
789  makeFaceCentresAndAreas(newPoints, fCtrs, fAreas, magfAreas);
790 
791  // Check cell volumes and calculate new cell centres
792  vectorField cellCtrs(nCells());
794 
795  makeCellCentresAndVols(fCtrs, fAreas, cellCtrs, cellVols);
796 
797  // Check cell volumes
798  bool error = checkCellVolumes
799  (
800  cellVols, // vols
801  report, // report
802  detailedReport, // detailedReport
803  nullptr // setPtr
804  );
805 
806 
807  // Check face areas
808  bool areaError = checkFaceAreas
809  (
810  fAreas,
811  report, // report
812  detailedReport, // detailedReport,
813  nullptr // setPtr
814  );
815  error = error || areaError;
816 
817 
818  // Check pyramid volumes
819  bool pyrVolError = checkFacePyramids
820  (
821  newPoints,
822  cellCtrs,
823  report, // report,
824  detailedReport, // detailedReport,
825  -small, // minPyrVol
826  nullptr // setPtr
827  );
828  error = error || pyrVolError;
829 
830 
831  // Check face non-orthogonality
832  bool nonOrthoError = checkFaceOrthogonality
833  (
834  fAreas,
835  cellCtrs,
836  report, // report
837  detailedReport, // detailedReport
838  nullptr // setPtr
839  );
840  error = error || nonOrthoError;
841 
842 
843  if (!error && (debug || report))
844  {
845  Pout<< "Mesh motion check OK." << endl;
846  }
847 
848  return error;
849 }
850 
851 
852 // ************************************************************************* //
bool checkFaceAreas(const vectorField &faceAreas, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative face areas.
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void makeFaceCentresAndAreas(const pointField &p, vectorField &fCtrs, vectorField &fAreas, scalarField &magfAreas) const
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1243
static PackedBoolList getInternalOrCoupledFaces(const polyMesh &)
Get per face whether it is internal or coupled.
Definition: syncTools.C:217
label nFaces() const
Unit conversion functions.
static tmp< scalarField > faceWeights(const polyMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate interpolation factors field.
uint8_t direction
Definition: direction.H:45
List< face > faceList
Definition: faceListFwd.H:43
label nCells() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
volVectorField vectorField(fieldObject, mesh)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
static PackedBoolList getMasterFaces(const polyMesh &)
Get per face whether it is uncoupled or a master of a.
Definition: syncTools.C:153
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
static PackedBoolList getInternalOrMasterFaces(const polyMesh &)
Get per face whether it is internal or a master of a.
Definition: syncTools.C:181
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:67
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:894
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
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 const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dimensionedScalar cos(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1237
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
static tmp< scalarField > faceSkewness(const polyMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
Definition: polyMeshTools.C:90
List< label > labelList
A List of labels.
Definition: labelList.H:56
void makeCellCentresAndVols(const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols) const
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
const vectorField & cellCentres() const
static tmp< scalarField > volRatio(const polyMesh &mesh, const scalarField &vol)
Generate volume ratio field.
bool checkCellVolumes(const scalarField &vols, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check for negative cell volumes.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross product operators.
Definition: Vector.H:57
const scalarField & cellVols
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const vectorField & faceCentres() const
static tmp< scalarField > faceOrthogonality(const polyMesh &mesh, const vectorField &fAreas, const vectorField &cellCtrs)
Generate orthogonality field. (1 for fully orthogonal, < 1 for.
Definition: polyMeshTools.C:35
#define WarningInFunction
Report a warning using Foam::Warning.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const vectorField & faceAreas() const
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
static scalar skewThreshold_
Skewness warning threshold.
static scalar nonOrthThreshold_
Non-orthogonality warning threshold in deg.
static tmp< scalarField > cellDeterminant(const primitiveMesh &mesh, const Vector< label > &directions, const vectorField &faceAreas, const PackedBoolList &internalOrCoupledFace)
Generate cell determinant field.
virtual bool checkMeshMotion(const pointField &newPoints, const bool report=false, const bool detailedReport=false) const
Check mesh motion for correctness given motion points.
const scalarField & cellVolumes() const
#define InfoInFunction
Report an information message using Foam::Info.