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