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 "polyMesh.H"
27 #include "polyMeshTools.H"
28 #include "unitConversion.H"
29 #include "syncTools.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
34 (
35  const bool report,
36  labelHashSet* setPtr
37 ) const
38 {
39  if (debug)
40  {
41  InfoInFunction << "Checking mesh non-orthogonality" << endl;
42  }
43 
44  const vectorField& fAreas = faceAreas();
45  const vectorField& cellCtrs = cellCentres();
46 
47  // Calculate orthogonality for all internal and coupled boundary faces
48  // (1 for uncoupled boundary faces)
50  (
51  *this,
52  fAreas,
53  cellCtrs
54  );
55  const scalarField& ortho = tortho.ref();
56 
57  // Severe nonorthogonality threshold
58  const scalar severeNonorthogonalityThreshold =
60 
61 
62  scalar minDDotS = great;
63  scalar sumDDotS = 0.0;
64  label nSummed = 0;
65  label severeNonOrth = 0;
66  label errorNonOrth = 0;
67 
68 
69  // Statistics only for internal and masters of coupled faces
71 
72  forAll(ortho, facei)
73  {
74  if (ortho[facei] < severeNonorthogonalityThreshold)
75  {
76  if (ortho[facei] > small)
77  {
78  if (setPtr)
79  {
80  setPtr->insert(facei);
81  }
82 
83  severeNonOrth++;
84  }
85  else
86  {
87  // Error : non-ortho too large
88  if (setPtr)
89  {
90  setPtr->insert(facei);
91  }
92 
93  errorNonOrth++;
94  }
95  }
96 
97  if (isMasterFace[facei])
98  {
99  minDDotS = min(minDDotS, ortho[facei]);
100  sumDDotS += ortho[facei];
101  nSummed++;
102  }
103  }
104 
105  reduce(minDDotS, minOp<scalar>());
106  reduce(sumDDotS, sumOp<scalar>());
107  reduce(nSummed, sumOp<label>());
108  reduce(severeNonOrth, sumOp<label>());
109  reduce(errorNonOrth, sumOp<label>());
110 
111  if (debug || report)
112  {
113  if (nSummed > 0)
114  {
115  if (debug || report)
116  {
117  Info<< " Mesh non-orthogonality Max: "
118  << radToDeg(::acos(min(1.0, max(-1.0, minDDotS))))
119  << " average: "
120  << radToDeg(::acos(min(1.0, max(-1.0, sumDDotS/nSummed))))
121  << endl;
122  }
123  }
124 
125  if (severeNonOrth > 0)
126  {
127  Info<< " *Number of severely non-orthogonal (> "
128  << primitiveMeshCheck::nonOrthThreshold << " degrees) faces: "
129  << severeNonOrth << "." << endl;
130  }
131  }
132 
133  if (errorNonOrth > 0)
134  {
135  if (debug || report)
136  {
137  Info<< " ***Number of non-orthogonality errors: "
138  << errorNonOrth << "." << endl;
139  }
140 
141  return true;
142  }
143  else
144  {
145  if (debug || report)
146  {
147  Info<< " Non-orthogonality check OK." << endl;
148  }
149 
150  return false;
151  }
152 }
153 
154 
156 (
157  const bool report,
158  labelHashSet* setPtr
159 ) const
160 {
161  if (debug)
162  {
163  InfoInFunction << "Checking face skewness" << endl;
164  }
165 
166  const pointField& points = this->points();
167  const vectorField& fCtrs = faceCentres();
168  const vectorField& fAreas = faceAreas();
169  const vectorField& cellCtrs = cellCentres();
170 
171  // Warn if the skew correction vector is more than skewWarning times
172  // larger than the face area vector
173 
175  (
176  *this,
177  points,
178  fCtrs,
179  fAreas,
180  cellCtrs
181  );
182  const scalarField& skew = tskew.ref();
183 
184  scalar maxSkew = max(skew);
185  label nWarnSkew = 0;
186 
187  // Statistics only for all faces except slave coupled faces
188  PackedBoolList isMasterFace(syncTools::getMasterFaces(*this));
189 
190  forAll(skew, facei)
191  {
192  // Check if the skewness vector is greater than the PN vector.
193  // This does not cause trouble but is a good indication of a poor mesh.
195  {
196  if (setPtr)
197  {
198  setPtr->insert(facei);
199  }
200 
201  if (isMasterFace[facei])
202  {
203  nWarnSkew++;
204  }
205  }
206  }
207 
208  reduce(maxSkew, maxOp<scalar>());
209  reduce(nWarnSkew, sumOp<label>());
210 
211  if (nWarnSkew > 0)
212  {
213  if (debug || report)
214  {
215  Info<< " ***Max skewness = " << maxSkew
216  << ", " << nWarnSkew << " highly skew faces detected"
217  " which may impair the quality of the results"
218  << endl;
219  }
220 
221  return true;
222  }
223  else
224  {
225  if (debug || report)
226  {
227  Info<< " Max skewness = " << maxSkew << " OK." << endl;
228  }
229 
230  return false;
231  }
232 }
233 
234 
236 (
237  const bool report,
238  const Vector<label>& directions,
239  labelHashSet* setPtr
240 ) const
241 {
242  // Check 1D/2Dness of edges. Gets passed the non-empty directions and
243  // checks all edges in the mesh whether they:
244  // - have no component in a non-empty direction or
245  // - are only in a single non-empty direction.
246  // Empty direction info is passed in as a vector of labels (synchronised)
247  // which are 1 if the direction is non-empty, 0 if it is.
248 
249  if (debug)
250  {
251  InfoInFunction << "Checking edge alignment" << endl;
252  }
253 
254  const pointField& p = points();
255 
256  label nDirs = 0;
257  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
258  {
259  if (directions[cmpt] == 1)
260  {
261  nDirs++;
262  }
263  else if (directions[cmpt] != 0)
264  {
266  << "directions should contain 0 or 1 but is now " << directions
267  << exit(FatalError);
268  }
269  }
270 
271  if (nDirs == vector::nComponents)
272  {
273  return false;
274  }
275 
276 
277  const faceList& fcs = faces();
278 
279  EdgeMap<label> edgesInError;
280 
281  forAll(fcs, facei)
282  {
283  const face& f = fcs[facei];
284 
285  forAll(f, fp)
286  {
287  label p0 = f[fp];
288  label p1 = f.nextLabel(fp);
289  if (p0 < p1)
290  {
291  vector d(p[p1]-p[p0]);
292  scalar magD = mag(d);
293 
294  if (magD > rootVSmall)
295  {
296  d /= magD;
297 
298  // Check how many empty directions are used by the edge.
299  label nEmptyDirs = 0;
300  label nNonEmptyDirs = 0;
301  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
302  {
303  if (mag(d[cmpt]) > 1e-6)
304  {
305  if (directions[cmpt] == 0)
306  {
307  nEmptyDirs++;
308  }
309  else
310  {
311  nNonEmptyDirs++;
312  }
313  }
314  }
315 
316  if (nEmptyDirs == 0)
317  {
318  // Purely in ok directions.
319  }
320  else if (nEmptyDirs == 1)
321  {
322  // Ok if purely in empty directions.
323  if (nNonEmptyDirs > 0)
324  {
325  edgesInError.insert(edge(p0, p1), facei);
326  }
327  }
328  else if (nEmptyDirs > 1)
329  {
330  // Always an error
331  edgesInError.insert(edge(p0, p1), facei);
332  }
333  }
334  }
335  }
336  }
337 
338  label nErrorEdges = returnReduce(edgesInError.size(), sumOp<label>());
339 
340  if (nErrorEdges > 0)
341  {
342  if (debug || report)
343  {
344  Info<< " ***Number of edges not aligned with or perpendicular to "
345  << "non-empty directions: " << nErrorEdges << endl;
346  }
347 
348  if (setPtr)
349  {
350  setPtr->resize(2*edgesInError.size());
351  forAllConstIter(EdgeMap<label>, edgesInError, iter)
352  {
353  setPtr->insert(iter.key()[0]);
354  setPtr->insert(iter.key()[1]);
355  }
356  }
357 
358  return true;
359  }
360  else
361  {
362  if (debug || report)
363  {
364  Info<< " All edges aligned with or perpendicular to "
365  << "non-empty directions." << endl;
366  }
367  return false;
368  }
369 }
370 
371 
373 (
374  const bool report,
375  labelHashSet* setPtr
376 ) const
377 {
378  const vectorField& faceAreas = this->faceAreas();
379  const Vector<label>& meshD = geometricD();
380 
381  const scalar warnDet = 1e-3;
382 
383  if (debug)
384  {
385  InfoInFunction << "Checking for under-determined cells" << endl;
386  }
387 
389  (
390  *this,
391  meshD,
392  faceAreas,
394  );
395  scalarField& cellDeterminant = tcellDeterminant.ref();
396 
397 
398  label nErrorCells = 0;
399  scalar minDet = min(cellDeterminant);
400  scalar sumDet = sum(cellDeterminant);
401 
402  forAll(cellDeterminant, celli)
403  {
404  if (cellDeterminant[celli] < warnDet)
405  {
406  if (setPtr)
407  {
408  setPtr->insert(celli);
409  }
410 
411  nErrorCells++;
412  }
413  }
414 
415  reduce(nErrorCells, sumOp<label>());
416  reduce(minDet, minOp<scalar>());
417  reduce(sumDet, sumOp<scalar>());
418  label nSummed = returnReduce(cellDeterminant.size(), sumOp<label>());
419 
420  if (debug || report)
421  {
422  if (nSummed > 0)
423  {
424  Info<< " Cell determinant (wellposedness) : minimum: " << minDet
425  << " average: " << sumDet/nSummed
426  << endl;
427  }
428  }
429 
430  if (nErrorCells > 0)
431  {
432  if (debug || report)
433  {
434  Info<< " ***Cells with small determinant (< "
435  << warnDet << ") found, number of cells: "
436  << nErrorCells << endl;
437  }
438 
439  return true;
440  }
441  else
442  {
443  if (debug || report)
444  {
445  Info<< " Cell determinant check OK." << endl;
446  }
447 
448  return false;
449  }
450 
451  return false;
452 }
453 
454 
456 (
457  const bool report,
458  const scalar minWeight,
459  labelHashSet* setPtr
460 ) const
461 {
462  if (debug)
463  {
464  InfoInFunction << "Checking for low face interpolation weights" << endl;
465  }
466 
467  const vectorField& fCtrs = faceCentres();
468  const vectorField& fAreas = faceAreas();
469  const vectorField& cellCtrs = this->cellCentres();
470 
472  (
473  *this,
474  fCtrs,
475  fAreas,
476  cellCtrs
477  );
478  scalarField& faceWght = tfaceWght.ref();
479 
480 
481  label nErrorFaces = 0;
482  scalar minDet = great;
483  scalar sumDet = 0.0;
484  label nSummed = 0;
485 
486  // Statistics only for internal and masters of coupled faces
488 
489  forAll(faceWght, facei)
490  {
491  if (faceWght[facei] < minWeight)
492  {
493  // Note: insert both sides of coupled faces
494  if (setPtr)
495  {
496  setPtr->insert(facei);
497  }
498 
499  nErrorFaces++;
500  }
501 
502  // Note: statistics only on master of coupled faces
503  if (isMasterFace[facei])
504  {
505  minDet = min(minDet, faceWght[facei]);
506  sumDet += faceWght[facei];
507  nSummed++;
508  }
509  }
510 
511  reduce(nErrorFaces, sumOp<label>());
512  reduce(minDet, minOp<scalar>());
513  reduce(sumDet, sumOp<scalar>());
514  reduce(nSummed, sumOp<label>());
515 
516  if (debug || report)
517  {
518  if (nSummed > 0)
519  {
520  Info<< " Face interpolation weight : minimum: " << minDet
521  << " average: " << sumDet/nSummed
522  << endl;
523  }
524  }
525 
526  if (nErrorFaces > 0)
527  {
528  if (debug || report)
529  {
530  Info<< " ***Faces with small interpolation weight (< " << minWeight
531  << ") found, number of faces: "
532  << nErrorFaces << endl;
533  }
534 
535  return true;
536  }
537  else
538  {
539  if (debug || report)
540  {
541  Info<< " Face interpolation weight check OK." << endl;
542  }
543 
544  return false;
545  }
546 
547  return false;
548 }
549 
550 
552 (
553  const bool report,
554  const scalar minRatio,
555  labelHashSet* setPtr
556 ) const
557 {
558  if (debug)
559  {
560  InfoInFunction << "Checking for volume ratio < " << minRatio << endl;
561  }
562 
563  const scalarField& cellVols = cellVolumes();
564 
566  scalarField& volRatio = tvolRatio.ref();
567 
568 
569  label nErrorFaces = 0;
570  scalar minDet = great;
571  scalar sumDet = 0.0;
572  label nSummed = 0;
573 
574  // Statistics only for internal and masters of coupled faces
576 
577  forAll(volRatio, facei)
578  {
579  if (volRatio[facei] < minRatio)
580  {
581  // Note: insert both sides of coupled faces
582  if (setPtr)
583  {
584  setPtr->insert(facei);
585  }
586 
587  nErrorFaces++;
588  }
589 
590  // Note: statistics only on master of coupled faces
591  if (isMasterFace[facei])
592  {
593  minDet = min(minDet, volRatio[facei]);
594  sumDet += volRatio[facei];
595  nSummed++;
596  }
597  }
598 
599  reduce(nErrorFaces, sumOp<label>());
600  reduce(minDet, minOp<scalar>());
601  reduce(sumDet, sumOp<scalar>());
602  reduce(nSummed, sumOp<label>());
603 
604  if (debug || report)
605  {
606  if (nSummed > 0)
607  {
608  Info<< " Face volume ratio : minimum: " << minDet
609  << " average: " << sumDet/nSummed
610  << endl;
611  }
612  }
613 
614  if (nErrorFaces > 0)
615  {
616  if (debug || report)
617  {
618  Info<< " ***Faces with small volume ratio (< " << minRatio
619  << ") found, number of faces: "
620  << nErrorFaces << endl;
621  }
622 
623  return true;
624  }
625  else
626  {
627  if (debug || report)
628  {
629  Info<< " Face volume ratio check OK." << endl;
630  }
631 
632  return false;
633  }
634 
635  return false;
636 }
637 
638 
639 // ************************************************************************* //
#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.
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:34
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.
const vectorField & cellCentres() const
const vectorField & faceAreas() const
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 skewThreshold
Skewness warning threshold.
scalar nonOrthThreshold
Non-orthogonality warning threshold in deg.
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
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)
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.