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