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-2023 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 "polyMeshCheck.H"
27 #include "polyMeshTools.H"
28 #include "unitConversion.H"
29 #include "syncTools.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
35 Foam::scalar Foam::polyMeshCheck::nonOrthThreshold = 70; // deg
38 
39 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
40 
42 (
43  const bool report,
44  labelHashSet* setPtr
45 ) const
46 {
47  if (debug)
48  {
49  InfoInFunction << "Checking mesh non-orthogonality" << endl;
50  }
51 
52  const vectorField& fAreas = faceAreas();
53  const vectorField& cellCtrs = cellCentres();
54 
55  // Calculate orthogonality for all internal and coupled boundary faces
56  // (1 for uncoupled boundary faces)
58  (
59  *this,
60  fAreas,
61  cellCtrs
62  );
63  const scalarField& ortho = tortho.ref();
64 
65  // Severe nonorthogonality threshold
66  const scalar severeNonorthogonalityThreshold =
68 
69 
70  scalar minDDotS = great;
71  scalar sumDDotS = 0.0;
72  label nSummed = 0;
73  label severeNonOrth = 0;
74  label errorNonOrth = 0;
75 
76 
77  // Statistics only for internal and masters of coupled faces
79 
80  forAll(ortho, facei)
81  {
82  if (ortho[facei] < severeNonorthogonalityThreshold)
83  {
84  if (ortho[facei] > small)
85  {
86  if (setPtr)
87  {
88  setPtr->insert(facei);
89  }
90 
91  severeNonOrth++;
92  }
93  else
94  {
95  // Error : non-ortho too large
96  if (setPtr)
97  {
98  setPtr->insert(facei);
99  }
100 
101  errorNonOrth++;
102  }
103  }
104 
105  if (isMasterFace[facei])
106  {
107  minDDotS = min(minDDotS, ortho[facei]);
108  sumDDotS += ortho[facei];
109  nSummed++;
110  }
111  }
112 
113  reduce(minDDotS, minOp<scalar>());
114  reduce(sumDDotS, sumOp<scalar>());
115  reduce(nSummed, sumOp<label>());
116  reduce(severeNonOrth, sumOp<label>());
117  reduce(errorNonOrth, sumOp<label>());
118 
119  if (debug || report)
120  {
121  if (nSummed > 0)
122  {
123  if (debug || report)
124  {
125  Info<< " Mesh non-orthogonality Max: "
126  << radToDeg(::acos(min(1.0, max(-1.0, minDDotS))))
127  << " average: "
128  << radToDeg(::acos(min(1.0, max(-1.0, sumDDotS/nSummed))))
129  << endl;
130  }
131  }
132 
133  if (severeNonOrth > 0)
134  {
135  Info<< " *Number of severely non-orthogonal (> "
136  << polyMeshCheck::nonOrthThreshold << " degrees) faces: "
137  << severeNonOrth << "." << endl;
138  }
139  }
140 
141  if (errorNonOrth > 0)
142  {
143  if (debug || report)
144  {
145  Info<< " ***Number of non-orthogonality errors: "
146  << errorNonOrth << "." << endl;
147  }
148 
149  return true;
150  }
151  else
152  {
153  if (debug || report)
154  {
155  Info<< " Non-orthogonality check OK." << endl;
156  }
157 
158  return false;
159  }
160 }
161 
162 
164 (
165  const bool report,
166  labelHashSet* setPtr
167 ) const
168 {
169  if (debug)
170  {
171  InfoInFunction << "Checking face skewness" << endl;
172  }
173 
174  const pointField& points = this->points();
175  const vectorField& fCtrs = faceCentres();
176  const vectorField& fAreas = faceAreas();
177  const vectorField& cellCtrs = cellCentres();
178 
179  // Warn if the skew correction vector is more than skewWarning times
180  // larger than the face area vector
181 
183  (
184  *this,
185  points,
186  fCtrs,
187  fAreas,
188  cellCtrs
189  );
190  const scalarField& skew = tskew.ref();
191 
192  scalar maxSkew = max(skew);
193  label nWarnSkew = 0;
194 
195  // Statistics only for all faces except slave coupled faces
196  PackedBoolList isMasterFace(syncTools::getMasterFaces(*this));
197 
198  forAll(skew, facei)
199  {
200  // Check if the skewness vector is greater than the PN vector.
201  // This does not cause trouble but is a good indication of a poor mesh.
202  if (skew[facei] > polyMeshCheck::skewThreshold)
203  {
204  if (setPtr)
205  {
206  setPtr->insert(facei);
207  }
208 
209  if (isMasterFace[facei])
210  {
211  nWarnSkew++;
212  }
213  }
214  }
215 
216  reduce(maxSkew, maxOp<scalar>());
217  reduce(nWarnSkew, sumOp<label>());
218 
219  if (nWarnSkew > 0)
220  {
221  if (debug || report)
222  {
223  Info<< " ***Max skewness = " << maxSkew
224  << ", " << nWarnSkew << " highly skew faces detected"
225  " which may impair the quality of the results"
226  << endl;
227  }
228 
229  return true;
230  }
231  else
232  {
233  if (debug || report)
234  {
235  Info<< " Max skewness = " << maxSkew << " OK." << endl;
236  }
237 
238  return false;
239  }
240 }
241 
242 
244 (
245  const bool report,
246  const Vector<label>& directions,
247  labelHashSet* setPtr
248 ) const
249 {
250  // Check 1D/2Dness of edges. Gets passed the non-empty directions and
251  // checks all edges in the mesh whether they:
252  // - have no component in a non-empty direction or
253  // - are only in a single non-empty direction.
254  // Empty direction info is passed in as a vector of labels (synchronised)
255  // which are 1 if the direction is non-empty, 0 if it is.
256 
257  if (debug)
258  {
259  InfoInFunction << "Checking edge alignment" << endl;
260  }
261 
262  const pointField& p = points();
263 
264  label nDirs = 0;
265  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
266  {
267  if (directions[cmpt] == 1)
268  {
269  nDirs++;
270  }
271  else if (directions[cmpt] != 0)
272  {
274  << "directions should contain 0 or 1 but is now " << directions
275  << exit(FatalError);
276  }
277  }
278 
279  if (nDirs == vector::nComponents)
280  {
281  return false;
282  }
283 
284 
285  const faceList& fcs = faces();
286 
287  EdgeMap<label> edgesInError;
288 
289  forAll(fcs, facei)
290  {
291  const face& f = fcs[facei];
292 
293  forAll(f, fp)
294  {
295  label p0 = f[fp];
296  label p1 = f.nextLabel(fp);
297  if (p0 < p1)
298  {
299  vector d(p[p1]-p[p0]);
300  scalar magD = mag(d);
301 
302  if (magD > rootVSmall)
303  {
304  d /= magD;
305 
306  // Check how many empty directions are used by the edge.
307  label nEmptyDirs = 0;
308  label nNonEmptyDirs = 0;
309  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
310  {
311  if (mag(d[cmpt]) > 1e-6)
312  {
313  if (directions[cmpt] == 0)
314  {
315  nEmptyDirs++;
316  }
317  else
318  {
319  nNonEmptyDirs++;
320  }
321  }
322  }
323 
324  if (nEmptyDirs == 0)
325  {
326  // Purely in ok directions.
327  }
328  else if (nEmptyDirs == 1)
329  {
330  // Ok if purely in empty directions.
331  if (nNonEmptyDirs > 0)
332  {
333  edgesInError.insert(edge(p0, p1), facei);
334  }
335  }
336  else if (nEmptyDirs > 1)
337  {
338  // Always an error
339  edgesInError.insert(edge(p0, p1), facei);
340  }
341  }
342  }
343  }
344  }
345 
346  label nErrorEdges = returnReduce(edgesInError.size(), sumOp<label>());
347 
348  if (nErrorEdges > 0)
349  {
350  if (debug || report)
351  {
352  Info<< " ***Number of edges not aligned with or perpendicular to "
353  << "non-empty directions: " << nErrorEdges << endl;
354  }
355 
356  if (setPtr)
357  {
358  setPtr->resize(2*edgesInError.size());
359  forAllConstIter(EdgeMap<label>, edgesInError, iter)
360  {
361  setPtr->insert(iter.key()[0]);
362  setPtr->insert(iter.key()[1]);
363  }
364  }
365 
366  return true;
367  }
368  else
369  {
370  if (debug || report)
371  {
372  Info<< " All edges aligned with or perpendicular to "
373  << "non-empty directions." << endl;
374  }
375  return false;
376  }
377 }
378 
379 
381 (
382  const bool report,
383  labelHashSet* setPtr
384 ) const
385 {
386  const vectorField& faceAreas = this->faceAreas();
387  const Vector<label>& meshD = geometricD();
388 
389  const scalar warnDet = 1e-3;
390 
391  if (debug)
392  {
393  InfoInFunction << "Checking for under-determined cells" << endl;
394  }
395 
397  (
398  *this,
399  meshD,
400  faceAreas,
402  );
403  scalarField& cellDeterminant = tcellDeterminant.ref();
404 
405 
406  label nErrorCells = 0;
407  scalar minDet = min(cellDeterminant);
408  scalar sumDet = sum(cellDeterminant);
409 
410  forAll(cellDeterminant, celli)
411  {
412  if (cellDeterminant[celli] < warnDet)
413  {
414  if (setPtr)
415  {
416  setPtr->insert(celli);
417  }
418 
419  nErrorCells++;
420  }
421  }
422 
423  reduce(nErrorCells, sumOp<label>());
424  reduce(minDet, minOp<scalar>());
425  reduce(sumDet, sumOp<scalar>());
426  label nSummed = returnReduce(cellDeterminant.size(), sumOp<label>());
427 
428  if (debug || report)
429  {
430  if (nSummed > 0)
431  {
432  Info<< " Cell determinant (wellposedness) : minimum: " << minDet
433  << " average: " << sumDet/nSummed
434  << endl;
435  }
436  }
437 
438  if (nErrorCells > 0)
439  {
440  if (debug || report)
441  {
442  Info<< " ***Cells with small determinant (< "
443  << warnDet << ") found, number of cells: "
444  << nErrorCells << endl;
445  }
446 
447  return true;
448  }
449  else
450  {
451  if (debug || report)
452  {
453  Info<< " Cell determinant check OK." << endl;
454  }
455 
456  return false;
457  }
458 
459  return false;
460 }
461 
462 
464 (
465  const bool report,
466  const scalar minWeight,
467  labelHashSet* setPtr
468 ) const
469 {
470  if (debug)
471  {
472  InfoInFunction << "Checking for low face interpolation weights" << endl;
473  }
474 
475  const vectorField& fCtrs = faceCentres();
476  const vectorField& fAreas = faceAreas();
477  const vectorField& cellCtrs = this->cellCentres();
478 
480  (
481  *this,
482  fCtrs,
483  fAreas,
484  cellCtrs
485  );
486  scalarField& faceWght = tfaceWght.ref();
487 
488 
489  label nErrorFaces = 0;
490  scalar minDet = great;
491  scalar sumDet = 0.0;
492  label nSummed = 0;
493 
494  // Statistics only for internal and masters of coupled faces
496 
497  forAll(faceWght, facei)
498  {
499  if (faceWght[facei] < minWeight)
500  {
501  // Note: insert both sides of coupled faces
502  if (setPtr)
503  {
504  setPtr->insert(facei);
505  }
506 
507  nErrorFaces++;
508  }
509 
510  // Note: statistics only on master of coupled faces
511  if (isMasterFace[facei])
512  {
513  minDet = min(minDet, faceWght[facei]);
514  sumDet += faceWght[facei];
515  nSummed++;
516  }
517  }
518 
519  reduce(nErrorFaces, sumOp<label>());
520  reduce(minDet, minOp<scalar>());
521  reduce(sumDet, sumOp<scalar>());
522  reduce(nSummed, sumOp<label>());
523 
524  if (debug || report)
525  {
526  if (nSummed > 0)
527  {
528  Info<< " Face interpolation weight : minimum: " << minDet
529  << " average: " << sumDet/nSummed
530  << endl;
531  }
532  }
533 
534  if (nErrorFaces > 0)
535  {
536  if (debug || report)
537  {
538  Info<< " ***Faces with small interpolation weight (< " << minWeight
539  << ") found, number of faces: "
540  << nErrorFaces << endl;
541  }
542 
543  return true;
544  }
545  else
546  {
547  if (debug || report)
548  {
549  Info<< " Face interpolation weight check OK." << endl;
550  }
551 
552  return false;
553  }
554 
555  return false;
556 }
557 
558 
560 (
561  const bool report,
562  const scalar minRatio,
563  labelHashSet* setPtr
564 ) const
565 {
566  if (debug)
567  {
568  InfoInFunction << "Checking for volume ratio < " << minRatio << endl;
569  }
570 
571  const scalarField& cellVols = cellVolumes();
572 
574  scalarField& volRatio = tvolRatio.ref();
575 
576 
577  label nErrorFaces = 0;
578  scalar minDet = great;
579  scalar sumDet = 0.0;
580  label nSummed = 0;
581 
582  // Statistics only for internal and masters of coupled faces
584 
585  forAll(volRatio, facei)
586  {
587  if (volRatio[facei] < minRatio)
588  {
589  // Note: insert both sides of coupled faces
590  if (setPtr)
591  {
592  setPtr->insert(facei);
593  }
594 
595  nErrorFaces++;
596  }
597 
598  // Note: statistics only on master of coupled faces
599  if (isMasterFace[facei])
600  {
601  minDet = min(minDet, volRatio[facei]);
602  sumDet += volRatio[facei];
603  nSummed++;
604  }
605  }
606 
607  reduce(nErrorFaces, sumOp<label>());
608  reduce(minDet, minOp<scalar>());
609  reduce(sumDet, sumOp<scalar>());
610  reduce(nSummed, sumOp<label>());
611 
612  if (debug || report)
613  {
614  if (nSummed > 0)
615  {
616  Info<< " Face volume ratio : minimum: " << minDet
617  << " average: " << sumDet/nSummed
618  << endl;
619  }
620  }
621 
622  if (nErrorFaces > 0)
623  {
624  if (debug || report)
625  {
626  Info<< " ***Faces with small volume ratio (< " << minRatio
627  << ") found, number of faces: "
628  << nErrorFaces << endl;
629  }
630 
631  return true;
632  }
633  else
634  {
635  if (debug || report)
636  {
637  Info<< " Face volume ratio check OK." << endl;
638  }
639 
640  return false;
641  }
642 
643  return false;
644 }
645 
646 
647 
648 bool Foam::polyMeshCheck::checkTopology(const polyMesh& mesh, const bool report)
649 {
650  label noFailedChecks = 0;
651 
652  if (mesh.checkPoints(report)) noFailedChecks++;
653  if (mesh.checkUpperTriangular(report)) noFailedChecks++;
654  if (mesh.checkCellsZipUp(report)) noFailedChecks++;
655  if (mesh.checkFaceVertices(report)) noFailedChecks++;
656  if (mesh.checkFaceFaces(report)) noFailedChecks++;
657 
658  if (noFailedChecks == 0)
659  {
660  if (report)
661  {
662  Info<< " Mesh topology OK." << endl;
663  }
664 
665  return false;
666  }
667  else
668  {
669  if (report)
670  {
671  Info<< " Failed " << noFailedChecks
672  << " mesh topology checks." << endl;
673  }
674 
675  return true;
676  }
677 }
678 
679 
680 bool Foam::polyMeshCheck::checkGeometry(const polyMesh& mesh, const bool report)
681 {
682  label noFailedChecks = 0;
683 
684  if (mesh.checkClosedBoundary(report)) noFailedChecks++;
685  if (mesh.checkClosedCells(report)) noFailedChecks++;
686  if (mesh.checkFaceAreas(report)) noFailedChecks++;
687  if (mesh.checkCellVolumes(report)) noFailedChecks++;
688  if (mesh.checkFaceOrthogonality(report)) noFailedChecks++;
689  if (mesh.checkFacePyramids(report)) noFailedChecks++;
690  if (mesh.checkFaceSkewness(report)) noFailedChecks++;
691 
692  if (noFailedChecks == 0)
693  {
694  if (report)
695  {
696  Info<< " Mesh geometry OK." << endl;
697  }
698 
699  return false;
700  }
701  else
702  {
703  if (report)
704  {
705  Info<< " Failed " << noFailedChecks
706  << " mesh geometry checks." << endl;
707  }
708 
709  return true;
710  }
711 }
712 
713 
714 bool Foam::polyMeshCheck::checkMesh(const polyMesh& mesh, const bool report)
715 {
716  const label noFailedChecks =
717  checkTopology(mesh, report)
718  + checkGeometry(mesh, report);
719 
720  if (noFailedChecks == 0)
721  {
722  if (report)
723  {
724  Info<< "Mesh OK." << endl;
725  }
726 
727  return false;
728  }
729  else
730  {
731  if (report)
732  {
733  Info<< " Failed " << noFailedChecks
734  << " mesh checks." << endl;
735  }
736 
737  return true;
738  }
739 }
740 
741 
742 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
void resize(const label newSize)
Resize the hash table for efficiency.
Definition: HashTable.C:432
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A bit-packed bool list.
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
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 point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
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 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
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.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual bool checkFaceSkewness(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face skewness.
virtual bool checkFaceOrthogonality(const bool report=false, labelHashSet *setPtr=nullptr) const
Check non-orthogonality.
Definition: polyMeshCheck.C:42
virtual bool checkEdgeAlignment(const bool report, const Vector< label > &directions, labelHashSet *setPtr) const
Check edge alignment for 1D/2D cases.
virtual bool checkFaceWeight(const bool report, const scalar minWeight=0.05, labelHashSet *setPtr=nullptr) const
Check for face weights.
virtual bool checkCellDeterminant(const bool report, labelHashSet *setPtr) const
virtual bool checkVolRatio(const bool report, const scalar minRatio=0.01, labelHashSet *setPtr=nullptr) const
Check for neighbouring cell volumes.
static tmp< scalarField > cellDeterminant(const primitiveMesh &mesh, const Vector< label > &directions, const vectorField &faceAreas, const PackedBoolList &internalOrCoupledFace)
Generate cell determinant field.
bool checkClosedCells(const bool report=false, labelHashSet *setPtr=nullptr, labelHashSet *highAspectSetPtr=nullptr, const Vector< label > &solutionD=Vector< label >::one) const
Check cells for closedness.
bool checkFacePyramids(const bool report=false, const scalar minPyrVol=-small, labelHashSet *setPtr=nullptr) const
Check face pyramid volume.
bool checkCellsZipUp(const bool report=false, labelHashSet *setPtr=nullptr) const
Check cell zip-up.
bool checkFaceAreas(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for negative face areas.
bool checkUpperTriangular(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face ordering.
const vectorField & cellCentres() const
bool checkClosedBoundary(const bool report=false) const
Check boundary for closedness.
bool checkCellVolumes(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for negative cell volumes.
bool checkFaceVertices(const bool report=false, labelHashSet *setPtr=nullptr) const
Check uniqueness of face vertices.
bool checkPoints(const bool report=false, labelHashSet *setPtr=nullptr) const
Check for unused points.
const vectorField & faceAreas() const
bool checkFaceFaces(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face-face connectivity.
static PackedBoolList getMasterFaces(const polyMesh &)
Get per face whether it is uncoupled or a master of a.
Definition: syncTools.C:153
static PackedBoolList getInternalOrCoupledFaces(const polyMesh &)
Get per face whether it is internal or coupled.
Definition: syncTools.C:217
static PackedBoolList getInternalOrMasterFaces(const polyMesh &)
Get per face whether it is internal or a master of a.
Definition: syncTools.C:181
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const pointField & points
#define InfoInFunction
Report an information message using Foam::Info.
scalar closedThreshold
Data to control mesh checking.
Definition: polyMeshCheck.C:33
bool checkTopology(const polyMesh &mesh, const bool report=false)
Check mesh topology for correctness.
scalar skewThreshold
Skewness warning threshold.
Definition: polyMeshCheck.C:36
bool checkMesh(const polyMesh &mesh, const bool report=false)
Check mesh for correctness. Returns false for no error.
scalar aspectThreshold
Aspect ratio warning threshold.
Definition: polyMeshCheck.C:34
scalar planarCosAngle
Threshold where faces are considered coplanar.
Definition: polyMeshCheck.C:37
scalar nonOrthThreshold
Non-orthogonality warning threshold in deg.
Definition: polyMeshCheck.C:35
bool checkGeometry(const polyMesh &mesh, const bool report=false)
Check mesh geometry (& implicitly topology) for correctness.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:105
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
messageStream Info
label checkGeometry(const polyMesh &mesh, const bool allGeometry, const autoPtr< surfaceWriter > &, const autoPtr< setWriter > &)
Check the geometry.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
label checkTopology(const polyMesh &mesh, const bool allTopology, const bool allGeometry, const autoPtr< surfaceWriter > &surfWriter, const autoPtr< Foam::setWriter > &setWriter)
error FatalError
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedTensor skew(const dimensionedTensor &dt)
uint8_t direction
Definition: direction.H:45
dimensionedScalar acos(const dimensionedScalar &ds)
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
labelList f(nPoints)
volScalarField & p
const scalarField & cellVols
Unit conversion functions.