polyMeshCheckQuality.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) 2011-2024 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"
28 #include "pyramidPointFaceRef.H"
29 #include "tetPointRef.H"
30 #include "syncTools.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace meshCheck
37 {
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
42 (
43  const primitiveMesh& mesh,
44  const bool report,
45  const scalar severeNonorthogonalityThreshold,
46  const label facei,
47  const vector& s, // face area vector
48  const vector& d, // cc-cc vector
49 
50  label& severeNonOrth,
51  label& errorNonOrth,
52  labelHashSet* setPtr
53 )
54 {
55  const scalar dDotS = (d & s)/(mag(d)*mag(s) + vSmall);
56 
57  if (dDotS < severeNonorthogonalityThreshold)
58  {
59  label nei = -1;
60 
61  if (mesh.isInternalFace(facei))
62  {
63  nei = mesh.faceNeighbour()[facei];
64  }
65 
66  if (dDotS > small)
67  {
68  if (report)
69  {
70  // Severe non-orthogonality but mesh still OK
71  Pout<< "Severe non-orthogonality for face " << facei
72  << " between cells " << mesh.faceOwner()[facei]
73  << " and " << nei
74  << ": Angle = "
75  << radToDeg(::acos(dDotS))
76  << " deg." << endl;
77  }
78 
79  severeNonOrth++;
80  }
81  else
82  {
83  // Non-orthogonality greater than 90 deg
84  if (report)
85  {
87  << "Severe non-orthogonality detected for face "
88  << facei
89  << " between cells " << mesh.faceOwner()[facei]
90  << " and " << nei
91  << ": Angle = "
92  << radToDeg(::acos(dDotS))
93  << " deg." << endl;
94  }
95 
96  errorNonOrth++;
97  }
98 
99  if (setPtr)
100  {
101  setPtr->insert(facei);
102  }
103  }
104 
105  return dDotS;
106 }
107 
108 
110 (
111  const primitiveMesh& mesh,
112  const bool report,
113  const scalar minTetQuality,
114  const pointField& p,
115  const label facei,
116  const point& fc, // face centre
117  const point& cc, // cell centre
118 
119  labelHashSet* setPtr
120 )
121 {
122  const face& f = mesh.faces()[facei];
123 
124  forAll(f, fp)
125  {
126  const scalar tetQual = tetPointRef
127  (
128  p[f[fp]],
129  p[f.nextLabel(fp)],
130  fc,
131  cc
132  ).quality();
133 
134  if (tetQual < minTetQuality)
135  {
136  if (report)
137  {
138  Pout<< "bool meshCheck::checkFaceTets("
139  << "const bool, const scalar, const pointField&"
140  << ", const pointField&"
141  << ", const labelList&, labelHashSet*) : "
142  << "face " << facei
143  << " has a triangle that points the wrong way."
144  << endl
145  << "Tet quality: " << tetQual
146  << " Face " << facei
147  << endl;
148  }
149 
150  if (setPtr)
151  {
152  setPtr->insert(facei);
153  }
154 
155  return true;
156  }
157  }
158 
159  return false;
160 }
161 
162 
164 (
165  const primitiveMesh& mesh,
166  const labelList& changedFaces
167 )
168 {
169  const labelList& own = mesh.faceOwner();
170  const labelList& nei = mesh.faceNeighbour();
171 
172  labelHashSet affectedCells(2*changedFaces.size());
173 
174  forAll(changedFaces, i)
175  {
176  const label facei = changedFaces[i];
177 
178  affectedCells.insert(own[facei]);
179 
180  if (mesh.isInternalFace(facei))
181  {
182  affectedCells.insert(nei[facei]);
183  }
184  }
185 
186  return affectedCells.toc();
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 } // End namespace meshCheck
193 } // End namespace Foam
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
198 (
199  const bool report,
200  const scalar orthWarn,
201  const polyMesh& mesh,
202  const vectorField& cellCentres,
203  const vectorField& faceAreas,
204  const labelList& checkFaces,
205  const List<labelPair>& baffles,
206  labelHashSet* setPtr
207 )
208 {
209  // for all internal and coupled faces check theat the d dot S product
210  // is positive
211 
212  const labelList& own = mesh.faceOwner();
213  const labelList& nei = mesh.faceNeighbour();
214  const polyBoundaryMesh& patches = mesh.boundaryMesh();
215 
216  // Severe nonorthogonality threshold
217  const scalar severeNonorthogonalityThreshold = ::cos(orthWarn);
218 
219  // Calculate coupled cell centre
220  pointField neiCc(mesh.nFaces() - mesh.nInternalFaces());
221 
222  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
223  {
224  neiCc[facei-mesh.nInternalFaces()] = cellCentres[own[facei]];
225  }
226 
227  syncTools::swapBoundaryFacePositions(mesh, neiCc);
228 
229  scalar minDDotS = great;
230 
231  scalar sumDDotS = 0;
232  label nDDotS = 0;
233 
234  label severeNonOrth = 0;
235 
236  label errorNonOrth = 0;
237 
238  forAll(checkFaces, i)
239  {
240  const label facei = checkFaces[i];
241 
242  const point& ownCc = cellCentres[own[facei]];
243 
244  if (mesh.isInternalFace(facei))
245  {
246  const scalar dDotS = checkNonOrtho
247  (
248  mesh,
249  report,
250  severeNonorthogonalityThreshold,
251  facei,
252  faceAreas[facei],
253  cellCentres[nei[facei]] - ownCc,
254 
255  severeNonOrth,
256  errorNonOrth,
257  setPtr
258  );
259 
260  if (dDotS < minDDotS)
261  {
262  minDDotS = dDotS;
263  }
264 
265  sumDDotS += dDotS;
266  nDDotS++;
267  }
268  else
269  {
270  const label patchi = patches.whichPatch(facei);
271 
272  if (patches[patchi].coupled())
273  {
274  const scalar dDotS = checkNonOrtho
275  (
276  mesh,
277  report,
278  severeNonorthogonalityThreshold,
279  facei,
280  faceAreas[facei],
281  neiCc[facei-mesh.nInternalFaces()] - ownCc,
282 
283  severeNonOrth,
284  errorNonOrth,
285  setPtr
286  );
287 
288  if (dDotS < minDDotS)
289  {
290  minDDotS = dDotS;
291  }
292 
293  sumDDotS += dDotS;
294  nDDotS++;
295  }
296  }
297  }
298 
299  forAll(baffles, i)
300  {
301  const label face0 = baffles[i].first();
302  const label face1 = baffles[i].second();
303 
304  const point& ownCc = cellCentres[own[face0]];
305 
306  const scalar dDotS = checkNonOrtho
307  (
308  mesh,
309  report,
310  severeNonorthogonalityThreshold,
311  face0,
312  faceAreas[face0],
313  cellCentres[own[face1]] - ownCc,
314 
315  severeNonOrth,
316  errorNonOrth,
317  setPtr
318  );
319 
320  if (dDotS < minDDotS)
321  {
322  minDDotS = dDotS;
323  }
324 
325  sumDDotS += dDotS;
326  nDDotS++;
327  }
328 
329  reduce(minDDotS, minOp<scalar>());
330  reduce(sumDDotS, sumOp<scalar>());
331  reduce(nDDotS, sumOp<label>());
332  reduce(severeNonOrth, sumOp<label>());
333  reduce(errorNonOrth, sumOp<label>());
334 
335  // Only report if there are some internal faces
336  if (nDDotS > 0)
337  {
338  if (report && minDDotS < severeNonorthogonalityThreshold)
339  {
340  Info<< "Number of non-orthogonality errors: " << errorNonOrth
341  << ". Number of severely non-orthogonal faces: "
342  << severeNonOrth << "." << endl;
343  }
344  }
345 
346  if (report)
347  {
348  if (nDDotS > 0)
349  {
350  Info<< "Mesh non-orthogonality Max: "
351  << radToDeg(::acos(minDDotS))
352  << " average: " << radToDeg(::acos(sumDDotS/nDDotS))
353  << endl;
354  }
355  }
356 
357  if (errorNonOrth > 0)
358  {
359  if (report)
360  {
362  << "Error in non-orthogonality detected" << endl;
363  }
364 
365  return true;
366  }
367  else
368  {
369  if (report)
370  {
371  Info<< "Non-orthogonality check OK.\n" << endl;
372  }
373 
374  return false;
375  }
376 }
377 
378 
380 (
381  const bool report,
382  const scalar minPyrVol,
383  const polyMesh& mesh,
384  const vectorField& cellCentres,
385  const pointField& p,
386  const labelList& checkFaces,
387  const List<labelPair>& baffles,
388  labelHashSet* setPtr
389 )
390 {
391  // check whether face area vector points to the cell with higher label
392  const labelList& own = mesh.faceOwner();
393  const labelList& nei = mesh.faceNeighbour();
394 
395  const faceList& f = mesh.faces();
396 
397  label nErrorPyrs = 0;
398 
399  forAll(checkFaces, i)
400  {
401  const label facei = checkFaces[i];
402 
403  // Create the owner pyramid - it will have negative volume
404  const scalar pyrVol = pyramidPointFaceRef
405  (
406  f[facei],
407  cellCentres[own[facei]]
408  ).mag(p);
409 
410  if (pyrVol > -minPyrVol)
411  {
412  if (report)
413  {
414  Pout<< "bool meshCheck::checkFacePyramids("
415  << "const bool, const scalar, const pointField&"
416  << ", const labelList&, labelHashSet*): "
417  << "face " << facei << " points the wrong way. " << endl
418  << "Pyramid volume: " << -pyrVol
419  << " Face " << f[facei] << " area: " << f[facei].mag(p)
420  << " Owner cell: " << own[facei] << endl
421  << "Owner cell vertex labels: "
422  << mesh.cells()[own[facei]].labels(f)
423  << endl;
424  }
425 
426 
427  if (setPtr)
428  {
429  setPtr->insert(facei);
430  }
431 
432  nErrorPyrs++;
433  }
434 
435  if (mesh.isInternalFace(facei))
436  {
437  // Create the neighbour pyramid - it will have positive volume
438  const scalar pyrVol =
439  pyramidPointFaceRef(f[facei], cellCentres[nei[facei]]).mag(p);
440 
441  if (pyrVol < minPyrVol)
442  {
443  if (report)
444  {
445  Pout<< "bool meshCheck::checkFacePyramids("
446  << "const bool, const scalar, const pointField&"
447  << ", const labelList&, labelHashSet*): "
448  << "face " << facei << " points the wrong way. " << endl
449  << "Pyramid volume: " << -pyrVol
450  << " Face " << f[facei] << " area: " << f[facei].mag(p)
451  << " Neighbour cell: " << nei[facei] << endl
452  << "Neighbour cell vertex labels: "
453  << mesh.cells()[nei[facei]].labels(f)
454  << endl;
455  }
456 
457  if (setPtr)
458  {
459  setPtr->insert(facei);
460  }
461 
462  nErrorPyrs++;
463  }
464  }
465  }
466 
467  forAll(baffles, i)
468  {
469  const label face0 = baffles[i].first();
470  const label face1 = baffles[i].second();
471 
472  const point& ownCc = cellCentres[own[face0]];
473 
474  // Create the owner pyramid - it will have negative volume
475  const scalar pyrVolOwn = pyramidPointFaceRef
476  (
477  f[face0],
478  ownCc
479  ).mag(p);
480 
481  if (pyrVolOwn > -minPyrVol)
482  {
483  if (report)
484  {
485  Pout<< "bool meshCheck::checkFacePyramids("
486  << "const bool, const scalar, const pointField&"
487  << ", const labelList&, labelHashSet*): "
488  << "face " << face0 << " points the wrong way. " << endl
489  << "Pyramid volume: " << -pyrVolOwn
490  << " Face " << f[face0] << " area: " << f[face0].mag(p)
491  << " Owner cell: " << own[face0] << endl
492  << "Owner cell vertex labels: "
493  << mesh.cells()[own[face0]].labels(f)
494  << endl;
495  }
496 
497 
498  if (setPtr)
499  {
500  setPtr->insert(face0);
501  }
502 
503  nErrorPyrs++;
504  }
505 
506  // Create the neighbour pyramid - it will have positive volume
507  const scalar pyrVolNbr =
508  pyramidPointFaceRef(f[face0], cellCentres[own[face1]]).mag(p);
509 
510  if (pyrVolNbr < minPyrVol)
511  {
512  if (report)
513  {
514  Pout<< "bool meshCheck::checkFacePyramids("
515  << "const bool, const scalar, const pointField&"
516  << ", const labelList&, labelHashSet*): "
517  << "face " << face0 << " points the wrong way. " << endl
518  << "Pyramid volume: " << -pyrVolNbr
519  << " Face " << f[face0] << " area: " << f[face0].mag(p)
520  << " Neighbour cell: " << own[face1] << endl
521  << "Neighbour cell vertex labels: "
522  << mesh.cells()[own[face1]].labels(f)
523  << endl;
524  }
525 
526  if (setPtr)
527  {
528  setPtr->insert(face0);
529  }
530 
531  nErrorPyrs++;
532  }
533  }
534 
535  reduce(nErrorPyrs, sumOp<label>());
536 
537  if (nErrorPyrs > 0)
538  {
539  if (report)
540  {
542  << "Error in face pyramids: faces pointing the wrong way."
543  << endl;
544  }
545 
546  return true;
547  }
548  else
549  {
550  if (report)
551  {
552  Info<< "Face pyramids OK.\n" << endl;
553  }
554 
555  return false;
556  }
557 }
558 
559 
561 (
562  const bool report,
563  const scalar minTetQuality,
564  const polyMesh& mesh,
565  const vectorField& cellCentres,
566  const vectorField& faceCentres,
567  const pointField& p,
568  const labelList& checkFaces,
569  const List<labelPair>& baffles,
570  labelHashSet* setPtr
571 )
572 {
573  // check whether decomposing each cell into tets results in
574  // positive volume, non-flat tets
575  const labelList& own = mesh.faceOwner();
576  const labelList& nei = mesh.faceNeighbour();
577  const polyBoundaryMesh& patches = mesh.boundaryMesh();
578 
579  // Calculate coupled cell centre
580  pointField neiCc(mesh.nFaces() - mesh.nInternalFaces());
581 
582  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
583  {
584  neiCc[facei - mesh.nInternalFaces()] = cellCentres[own[facei]];
585  }
586 
587  syncTools::swapBoundaryFacePositions(mesh, neiCc);
588 
589  label nErrorTets = 0;
590 
591  forAll(checkFaces, i)
592  {
593  const label facei = checkFaces[i];
594 
595  // Create the owner pyramid - note: exchange cell and face centre
596  // to get positive volume.
597  bool tetError = checkFaceTet
598  (
599  mesh,
600  report,
601  minTetQuality,
602  p,
603  facei,
604  cellCentres[own[facei]], // face centre
605  faceCentres[facei], // cell centre
606  setPtr
607  );
608 
609  if (tetError)
610  {
611  nErrorTets++;
612  }
613 
614  if (mesh.isInternalFace(facei))
615  {
616  // Create the neighbour tets - they will have positive volume
617  bool tetError = checkFaceTet
618  (
619  mesh,
620  report,
621  minTetQuality,
622  p,
623  facei,
624  faceCentres[facei], // face centre
625  cellCentres[nei[facei]], // cell centre
626  setPtr
627  );
628 
629  if (tetError)
630  {
631  nErrorTets++;
632  }
633 
634  if
635  (
636  polyMeshTetDecomposition::findSharedBasePoint
637  (
638  mesh,
639  facei,
640  minTetQuality,
641  report
642  ) == -1
643  )
644  {
645  if (setPtr)
646  {
647  setPtr->insert(facei);
648  }
649 
650  nErrorTets++;
651  }
652  }
653  else
654  {
655  const label patchi = patches.whichPatch(facei);
656 
657  if (patches[patchi].coupled())
658  {
659  if
660  (
661  polyMeshTetDecomposition::findSharedBasePoint
662  (
663  mesh,
664  facei,
665  neiCc[facei - mesh.nInternalFaces()],
666  minTetQuality,
667  report
668  ) == -1
669  )
670  {
671  if (setPtr)
672  {
673  setPtr->insert(facei);
674  }
675 
676  nErrorTets++;
677  }
678  }
679  else
680  {
681  if
682  (
683  polyMeshTetDecomposition::findBasePoint
684  (
685  mesh,
686  facei,
687  minTetQuality,
688  report
689  ) == -1
690  )
691  {
692  if (setPtr)
693  {
694  setPtr->insert(facei);
695  }
696 
697  nErrorTets++;
698  }
699  }
700  }
701  }
702 
703  forAll(baffles, i)
704  {
705  const label face0 = baffles[i].first();
706  const label face1 = baffles[i].second();
707 
708  bool tetError = checkFaceTet
709  (
710  mesh,
711  report,
712  minTetQuality,
713  p,
714  face0,
715  cellCentres[own[face0]], // face centre
716  faceCentres[face0], // cell centre
717  setPtr
718  );
719 
720  if (tetError)
721  {
722  nErrorTets++;
723  }
724 
725  // Create the neighbour tets - they will have positive volume
726  tetError = checkFaceTet
727  (
728  mesh,
729  report,
730  minTetQuality,
731  p,
732  face0,
733  faceCentres[face0], // face centre
734  cellCentres[own[face1]], // cell centre
735  setPtr
736  );
737 
738  if (tetError)
739  {
740  nErrorTets++;
741  }
742 
743  if
744  (
745  polyMeshTetDecomposition::findSharedBasePoint
746  (
747  mesh,
748  face0,
749  cellCentres[own[face1]],
750  minTetQuality,
751  report
752  ) == -1
753  )
754  {
755  if (setPtr)
756  {
757  setPtr->insert(face0);
758  }
759 
760  nErrorTets++;
761  }
762  }
763 
764  reduce(nErrorTets, sumOp<label>());
765 
766  if (nErrorTets > 0)
767  {
768  if (report)
769  {
771  << "Error in face decomposition: negative tets."
772  << endl;
773  }
774 
775  return true;
776  }
777  else
778  {
779  if (report)
780  {
781  Info<< "Face tets OK.\n" << endl;
782  }
783 
784  return false;
785  }
786 }
787 
788 
790 (
791  const bool report,
792  const scalar internalSkew,
793  const scalar boundarySkew,
794  const polyMesh& mesh,
795  const pointField& points,
796  const vectorField& cellCentres,
797  const vectorField& faceCentres,
798  const vectorField& faceAreas,
799  const labelList& checkFaces,
800  const List<labelPair>& baffles,
801  labelHashSet* setPtr
802 )
803 {
804  // Warn if the skew correction vector is more than skew times
805  // larger than the face area vector
806 
807  const labelList& own = mesh.faceOwner();
808  const labelList& nei = mesh.faceNeighbour();
809  const polyBoundaryMesh& patches = mesh.boundaryMesh();
810 
811  // Calculate coupled cell centre
812  pointField neiCc;
813  syncTools::swapBoundaryCellPositions(mesh, cellCentres, neiCc);
814 
815  scalar maxSkew = 0;
816 
817  label nWarnSkew = 0;
818 
819  forAll(checkFaces, i)
820  {
821  const label facei = checkFaces[i];
822 
823  if (mesh.isInternalFace(facei))
824  {
825  const scalar skewness = meshCheck::faceSkewness
826  (
827  mesh,
828  points,
829  faceCentres,
830  faceAreas,
831 
832  facei,
833  cellCentres[own[facei]],
834  cellCentres[nei[facei]]
835  );
836 
837  // Check if the skewness vector is greater than the PN vector.
838  // This does not cause trouble but is a good indication of a poor
839  // mesh.
840  if (skewness > internalSkew)
841  {
842  if (report)
843  {
844  Pout<< "Severe skewness for face " << facei
845  << " skewness = " << skewness << endl;
846  }
847 
848  if (setPtr)
849  {
850  setPtr->insert(facei);
851  }
852 
853  nWarnSkew++;
854  }
855 
856  maxSkew = max(maxSkew, skewness);
857  }
858  else if (patches[patches.whichPatch(facei)].coupled())
859  {
860  const scalar skewness = meshCheck::faceSkewness
861  (
862  mesh,
863  points,
864  faceCentres,
865  faceAreas,
866 
867  facei,
868  cellCentres[own[facei]],
869  neiCc[facei-mesh.nInternalFaces()]
870  );
871 
872  // Check if the skewness vector is greater than the PN vector.
873  // This does not cause trouble but is a good indication of a poor
874  // mesh.
875  if (skewness > internalSkew)
876  {
877  if (report)
878  {
879  Pout<< "Severe skewness for coupled face " << facei
880  << " skewness = " << skewness << endl;
881  }
882 
883  if (setPtr)
884  {
885  setPtr->insert(facei);
886  }
887 
888  nWarnSkew++;
889  }
890 
891  maxSkew = max(maxSkew, skewness);
892  }
893  else
894  {
895  const scalar skewness = meshCheck::boundaryFaceSkewness
896  (
897  mesh,
898  points,
899  faceCentres,
900  faceAreas,
901 
902  facei,
903  cellCentres[own[facei]]
904  );
905 
906 
907  // Check if the skewness vector is greater than the PN vector.
908  // This does not cause trouble but is a good indication of a poor
909  // mesh.
910  if (skewness > boundarySkew)
911  {
912  if (report)
913  {
914  Pout<< "Severe skewness for boundary face " << facei
915  << " skewness = " << skewness << endl;
916  }
917 
918  if (setPtr)
919  {
920  setPtr->insert(facei);
921  }
922 
923  nWarnSkew++;
924  }
925 
926  maxSkew = max(maxSkew, skewness);
927  }
928  }
929 
930  forAll(baffles, i)
931  {
932  const label face0 = baffles[i].first();
933  const label face1 = baffles[i].second();
934 
935  const point& ownCc = cellCentres[own[face0]];
936  const point& neiCc = cellCentres[own[face1]];
937 
938  const scalar skewness = meshCheck::faceSkewness
939  (
940  mesh,
941  points,
942  faceCentres,
943  faceAreas,
944 
945  face0,
946  ownCc,
947  neiCc
948  );
949 
950  // Check if the skewness vector is greater than the PN vector.
951  // This does not cause trouble but is a good indication of a poor
952  // mesh.
953  if (skewness > internalSkew)
954  {
955  if (report)
956  {
957  Pout<< "Severe skewness for face " << face0
958  << " skewness = " << skewness << endl;
959  }
960 
961  if (setPtr)
962  {
963  setPtr->insert(face0);
964  }
965 
966  nWarnSkew++;
967  }
968 
969  maxSkew = max(maxSkew, skewness);
970  }
971 
972 
973  reduce(maxSkew, maxOp<scalar>());
974  reduce(nWarnSkew, sumOp<label>());
975 
976  if (nWarnSkew > 0)
977  {
978  if (report)
979  {
981  << 100*maxSkew
982  << " percent.\nThis may impair the quality of the result." << nl
983  << nWarnSkew << " highly skew faces detected."
984  << endl;
985  }
986 
987  return true;
988  }
989  else
990  {
991  if (report)
992  {
993  Info<< "Max skewness = " << 100*maxSkew
994  << " percent. Face skewness OK.\n" << endl;
995  }
996 
997  return false;
998  }
999 }
1000 
1001 
1003 (
1004  const bool report,
1005  const scalar warnWeight,
1006  const polyMesh& mesh,
1007  const vectorField& cellCentres,
1008  const vectorField& faceCentres,
1009  const vectorField& faceAreas,
1010  const labelList& checkFaces,
1011  const List<labelPair>& baffles,
1012  labelHashSet* setPtr
1013 )
1014 {
1015  // Warn if the delta factor (0..1) is too large.
1016 
1017  const labelList& own = mesh.faceOwner();
1018  const labelList& nei = mesh.faceNeighbour();
1019  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1020 
1021  // Calculate coupled cell centre
1022  pointField neiCc(mesh.nFaces()-mesh.nInternalFaces());
1023 
1024  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
1025  {
1026  neiCc[facei-mesh.nInternalFaces()] = cellCentres[own[facei]];
1027  }
1028  syncTools::swapBoundaryFacePositions(mesh, neiCc);
1029 
1030 
1031  scalar minWeight = great;
1032 
1033  label nWarnWeight = 0;
1034 
1035  forAll(checkFaces, i)
1036  {
1037  const label facei = checkFaces[i];
1038 
1039  const point& fc = faceCentres[facei];
1040  const vector& fa = faceAreas[facei];
1041 
1042  const scalar dOwn = mag(fa & (fc-cellCentres[own[facei]]));
1043 
1044  if (mesh.isInternalFace(facei))
1045  {
1046  const scalar dNei = mag(fa & (cellCentres[nei[facei]]-fc));
1047  const scalar weight = min(dNei, dOwn)/(dNei + dOwn + vSmall);
1048 
1049  if (weight < warnWeight)
1050  {
1051  if (report)
1052  {
1053  Pout<< "Small weighting factor for face " << facei
1054  << " weight = " << weight << endl;
1055  }
1056 
1057  if (setPtr)
1058  {
1059  setPtr->insert(facei);
1060  }
1061 
1062  nWarnWeight++;
1063  }
1064 
1065  minWeight = min(minWeight, weight);
1066  }
1067  else
1068  {
1069  const label patchi = patches.whichPatch(facei);
1070 
1071  if (patches[patchi].coupled())
1072  {
1073  const scalar dNei =
1074  mag(fa & (neiCc[facei-mesh.nInternalFaces()]-fc));
1075  const scalar weight = min(dNei, dOwn)/(dNei + dOwn + vSmall);
1076 
1077  if (weight < warnWeight)
1078  {
1079  if (report)
1080  {
1081  Pout<< "Small weighting factor for face " << facei
1082  << " weight = " << weight << endl;
1083  }
1084 
1085  if (setPtr)
1086  {
1087  setPtr->insert(facei);
1088  }
1089 
1090  nWarnWeight++;
1091  }
1092 
1093  minWeight = min(minWeight, weight);
1094  }
1095  }
1096  }
1097 
1098  forAll(baffles, i)
1099  {
1100  const label face0 = baffles[i].first();
1101  const label face1 = baffles[i].second();
1102 
1103  const point& ownCc = cellCentres[own[face0]];
1104  const point& fc = faceCentres[face0];
1105  const vector& fa = faceAreas[face0];
1106 
1107  const scalar dOwn = mag(fa & (fc-ownCc));
1108  const scalar dNei = mag(fa & (cellCentres[own[face1]]-fc));
1109  const scalar weight = min(dNei, dOwn)/(dNei + dOwn + vSmall);
1110 
1111  if (weight < warnWeight)
1112  {
1113  if (report)
1114  {
1115  Pout<< "Small weighting factor for face " << face0
1116  << " weight = " << weight << endl;
1117  }
1118 
1119  if (setPtr)
1120  {
1121  setPtr->insert(face0);
1122  }
1123 
1124  nWarnWeight++;
1125  }
1126 
1127  minWeight = min(minWeight, weight);
1128  }
1129 
1130  reduce(minWeight, minOp<scalar>());
1131  reduce(nWarnWeight, sumOp<label>());
1132 
1133  if (minWeight < warnWeight)
1134  {
1135  if (report)
1136  {
1138  << minWeight << '.' << nl
1139  << nWarnWeight << " faces with small weights detected."
1140  << endl;
1141  }
1142 
1143  return true;
1144  }
1145  else
1146  {
1147  if (report)
1148  {
1149  Info<< "Min weight = " << minWeight
1150  << ". Weights OK.\n" << endl;
1151  }
1152 
1153  return false;
1154  }
1155 }
1156 
1157 
1159 (
1160  const bool report,
1161  const scalar warnRatio,
1162  const polyMesh& mesh,
1163  const scalarField& cellVolumes,
1164  const labelList& checkFaces,
1165  const List<labelPair>& baffles,
1166  labelHashSet* setPtr
1167 )
1168 {
1169  // Warn if the volume ratio between neighbouring cells is too large
1170 
1171  const labelList& own = mesh.faceOwner();
1172  const labelList& nei = mesh.faceNeighbour();
1173  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1174 
1175  // Calculate coupled cell vol
1176  scalarField neiVols(mesh.nFaces()-mesh.nInternalFaces());
1177 
1178  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
1179  {
1180  neiVols[facei-mesh.nInternalFaces()] = cellVolumes[own[facei]];
1181  }
1182  syncTools::swapBoundaryFaceList(mesh, neiVols);
1183 
1184 
1185  scalar minRatio = great;
1186 
1187  label nWarnRatio = 0;
1188 
1189  forAll(checkFaces, i)
1190  {
1191  const label facei = checkFaces[i];
1192 
1193  const scalar ownVol = mag(cellVolumes[own[facei]]);
1194 
1195  scalar neiVol = -great;
1196 
1197  if (mesh.isInternalFace(facei))
1198  {
1199  neiVol = mag(cellVolumes[nei[facei]]);
1200  }
1201  else
1202  {
1203  const label patchi = patches.whichPatch(facei);
1204 
1205  if (patches[patchi].coupled())
1206  {
1207  neiVol = mag(neiVols[facei-mesh.nInternalFaces()]);
1208  }
1209  }
1210 
1211  if (neiVol >= 0)
1212  {
1213  const scalar ratio =
1214  min(ownVol, neiVol)/(max(ownVol, neiVol) + vSmall);
1215 
1216  if (ratio < warnRatio)
1217  {
1218  if (report)
1219  {
1220  Pout<< "Small ratio for face " << facei
1221  << " ratio = " << ratio << endl;
1222  }
1223 
1224  if (setPtr)
1225  {
1226  setPtr->insert(facei);
1227  }
1228 
1229  nWarnRatio++;
1230  }
1231 
1232  minRatio = min(minRatio, ratio);
1233  }
1234  }
1235 
1236  forAll(baffles, i)
1237  {
1238  const label face0 = baffles[i].first();
1239  const label face1 = baffles[i].second();
1240 
1241  const scalar ownVol = mag(cellVolumes[own[face0]]);
1242 
1243  const scalar neiVol = mag(cellVolumes[own[face1]]);
1244 
1245  if (neiVol >= 0)
1246  {
1247  const scalar ratio =
1248  min(ownVol, neiVol)/(max(ownVol, neiVol) + vSmall);
1249 
1250  if (ratio < warnRatio)
1251  {
1252  if (report)
1253  {
1254  Pout<< "Small ratio for face " << face0
1255  << " ratio = " << ratio << endl;
1256  }
1257 
1258  if (setPtr)
1259  {
1260  setPtr->insert(face0);
1261  }
1262 
1263  nWarnRatio++;
1264  }
1265 
1266  minRatio = min(minRatio, ratio);
1267  }
1268  }
1269 
1270  reduce(minRatio, minOp<scalar>());
1271  reduce(nWarnRatio, sumOp<label>());
1272 
1273  if (minRatio < warnRatio)
1274  {
1275  if (report)
1276  {
1278  << minRatio << '.' << nl
1279  << nWarnRatio << " faces with small ratios detected."
1280  << endl;
1281  }
1282 
1283  return true;
1284  }
1285  else
1286  {
1287  if (report)
1288  {
1289  Info<< "Min ratio = " << minRatio
1290  << ". Ratios OK.\n" << endl;
1291  }
1292 
1293  return false;
1294  }
1295 }
1296 
1297 
1299 (
1300  const bool report,
1301  const scalar maxConcave,
1302  const polyMesh& mesh,
1303  const vectorField& faceAreas,
1304  const pointField& p,
1305  const labelList& checkFaces,
1306  labelHashSet* setPtr
1307 )
1308 {
1309  if (maxConcave < -small || maxConcave > degToRad(180)+small)
1310  {
1312  << "maxConcave should be [0..180] degrees but is "
1313  << radToDeg(maxConcave) << abort(FatalError);
1314  }
1315 
1316  const scalar maxSin = Foam::sin(maxConcave);
1317 
1318  const faceList& fcs = mesh.faces();
1319 
1320  scalar maxEdgeSin = 0.0;
1321 
1322  label nConcave = 0;
1323 
1324  label errorFacei = -1;
1325 
1326  forAll(checkFaces, i)
1327  {
1328  const label facei = checkFaces[i];
1329 
1330  const face& f = fcs[facei];
1331 
1332  vector faceNormal = faceAreas[facei];
1333  faceNormal /= mag(faceNormal) + vSmall;
1334 
1335  // Get edge from f[0] to f[size-1];
1336  vector ePrev(p[f.first()] - p[f.last()]);
1337  scalar magEPrev = mag(ePrev);
1338  ePrev /= magEPrev + vSmall;
1339 
1340  forAll(f, fp0)
1341  {
1342  // Get vertex after fp
1343  const label fp1 = f.fcIndex(fp0);
1344 
1345  // Normalised vector between two consecutive points
1346  vector e10(p[f[fp1]] - p[f[fp0]]);
1347  const scalar magE10 = mag(e10);
1348  e10 /= magE10 + vSmall;
1349 
1350  if (magEPrev > small && magE10 > small)
1351  {
1352  vector edgeNormal = ePrev ^ e10;
1353  const scalar magEdgeNormal = mag(edgeNormal);
1354 
1355  if (magEdgeNormal < maxSin)
1356  {
1357  // Edges (almost) aligned -> face is ok.
1358  }
1359  else
1360  {
1361  // Check normal
1362  edgeNormal /= magEdgeNormal;
1363 
1364  if ((edgeNormal & faceNormal) < small)
1365  {
1366  if (facei != errorFacei)
1367  {
1368  // Count only one error per face.
1369  errorFacei = facei;
1370  nConcave++;
1371  }
1372 
1373  if (setPtr)
1374  {
1375  setPtr->insert(facei);
1376  }
1377 
1378  maxEdgeSin = max(maxEdgeSin, magEdgeNormal);
1379  }
1380  }
1381  }
1382 
1383  ePrev = e10;
1384  magEPrev = magE10;
1385  }
1386  }
1387 
1388  reduce(nConcave, sumOp<label>());
1389  reduce(maxEdgeSin, maxOp<scalar>());
1390 
1391  if (report)
1392  {
1393  if (maxEdgeSin > small)
1394  {
1395  const scalar maxConcaveDegr =
1396  radToDeg(Foam::asin(Foam::min(1.0, maxEdgeSin)));
1397 
1398  Info<< "There are " << nConcave
1399  << " faces with concave angles between consecutive"
1400  << " edges. Max concave angle = " << maxConcaveDegr
1401  << " degrees.\n" << endl;
1402  }
1403  else
1404  {
1405  Info<< "All angles in faces are convex or less than "
1406  << radToDeg(maxConcave) << " degrees concave.\n" << endl;
1407  }
1408  }
1409 
1410  if (nConcave > 0)
1411  {
1412  if (report)
1413  {
1415  << nConcave << " face points with severe concave angle (> "
1416  << radToDeg(maxConcave) << " deg) found.\n"
1417  << endl;
1418  }
1419 
1420  return true;
1421  }
1422  else
1423  {
1424  return false;
1425  }
1426 }
1427 
1428 
1430 (
1431  const bool report,
1432  const scalar minTwist,
1433  const polyMesh& mesh,
1434  const vectorField& cellCentres,
1435  const vectorField& faceAreas,
1436  const vectorField& faceCentres,
1437  const pointField& p,
1438  const labelList& checkFaces,
1439  labelHashSet* setPtr
1440 )
1441 {
1442  if (minTwist < -1-small || minTwist > 1+small)
1443  {
1445  << "minTwist should be [-1..1] but is now " << minTwist
1446  << abort(FatalError);
1447  }
1448 
1449 
1450  const faceList& fcs = mesh.faces();
1451 
1452  label nWarped = 0;
1453 
1454  const labelList& own = mesh.faceOwner();
1455  const labelList& nei = mesh.faceNeighbour();
1456  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1457 
1458  // Calculate coupled cell centre
1459  pointField neiCc(mesh.nFaces()-mesh.nInternalFaces());
1460 
1461  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
1462  {
1463  neiCc[facei-mesh.nInternalFaces()] = cellCentres[own[facei]];
1464  }
1465  syncTools::swapBoundaryFacePositions(mesh, neiCc);
1466 
1467  forAll(checkFaces, i)
1468  {
1469  const label facei = checkFaces[i];
1470 
1471  const face& f = fcs[facei];
1472 
1473  if (f.size() > 3)
1474  {
1475  vector nf(Zero);
1476 
1477  if (mesh.isInternalFace(facei))
1478  {
1479  nf = cellCentres[nei[facei]] - cellCentres[own[facei]];
1480  nf /= mag(nf) + vSmall;
1481  }
1482  else if (patches[patches.whichPatch(facei)].coupled())
1483  {
1484  nf =
1485  neiCc[facei-mesh.nInternalFaces()]
1486  - cellCentres[own[facei]];
1487  nf /= mag(nf) + vSmall;
1488  }
1489  else
1490  {
1491  nf = faceCentres[facei] - cellCentres[own[facei]];
1492  nf /= mag(nf) + vSmall;
1493  }
1494 
1495  if (nf != vector::zero)
1496  {
1497  const point& fc = faceCentres[facei];
1498 
1499  forAll(f, fpI)
1500  {
1501  vector triArea
1502  (
1503  triPointRef
1504  (
1505  p[f[fpI]],
1506  p[f.nextLabel(fpI)],
1507  fc
1508  ).area()
1509  );
1510 
1511  const scalar magTri = mag(triArea);
1512 
1513  if (magTri > vSmall && ((nf & triArea/magTri) < minTwist))
1514  {
1515  nWarped++;
1516 
1517  if (setPtr)
1518  {
1519  setPtr->insert(facei);
1520  }
1521 
1522  break;
1523  }
1524  }
1525  }
1526  }
1527  }
1528 
1529  reduce(nWarped, sumOp<label>());
1530 
1531  if (report)
1532  {
1533  if (nWarped> 0)
1534  {
1535  Info<< "There are " << nWarped
1536  << " faces with cosine of the angle"
1537  << " between triangle normal and face normal less than "
1538  << minTwist << nl << endl;
1539  }
1540  else
1541  {
1542  Info<< "All faces are flat in that the cosine of the angle"
1543  << " between triangle normal and face normal less than "
1544  << minTwist << nl << endl;
1545  }
1546  }
1547 
1548  if (nWarped > 0)
1549  {
1550  if (report)
1551  {
1553  << nWarped << " faces with severe warpage "
1554  << "(cosine of the angle between triangle normal and "
1555  << "face normal < " << minTwist << ") found.\n"
1556  << endl;
1557  }
1558 
1559  return true;
1560  }
1561  else
1562  {
1563  return false;
1564  }
1565 }
1566 
1567 
1569 (
1570  const bool report,
1571  const scalar minTwist,
1572  const polyMesh& mesh,
1573  const vectorField& faceAreas,
1574  const vectorField& faceCentres,
1575  const pointField& p,
1576  const labelList& checkFaces,
1577  labelHashSet* setPtr
1578 )
1579 {
1580  if (minTwist < -1-small || minTwist > 1+small)
1581  {
1583  << "minTwist should be [-1..1] but is now " << minTwist
1584  << abort(FatalError);
1585  }
1586 
1587  const faceList& fcs = mesh.faces();
1588 
1589  label nWarped = 0;
1590 
1591  forAll(checkFaces, i)
1592  {
1593  const label facei = checkFaces[i];
1594 
1595  const face& f = fcs[facei];
1596 
1597  if (f.size() > 3)
1598  {
1599  const point& fc = faceCentres[facei];
1600 
1601  // Find starting triangle (at startFp) with non-zero area
1602  label startFp = -1;
1603  vector prevN;
1604 
1605  forAll(f, fp)
1606  {
1607  prevN = triPointRef
1608  (
1609  p[f[fp]],
1610  p[f.nextLabel(fp)],
1611  fc
1612  ).area();
1613 
1614  const scalar magTri = mag(prevN);
1615 
1616  if (magTri > vSmall)
1617  {
1618  startFp = fp;
1619  prevN /= magTri;
1620  break;
1621  }
1622  }
1623 
1624  if (startFp != -1)
1625  {
1626  label fp = startFp;
1627 
1628  do
1629  {
1630  fp = f.fcIndex(fp);
1631 
1632  vector triN
1633  (
1634  triPointRef
1635  (
1636  p[f[fp]],
1637  p[f.nextLabel(fp)],
1638  fc
1639  ).area()
1640  );
1641  const scalar magTri = mag(triN);
1642 
1643  if (magTri > vSmall)
1644  {
1645  triN /= magTri;
1646 
1647  if ((prevN & triN) < minTwist)
1648  {
1649  nWarped++;
1650 
1651  if (setPtr)
1652  {
1653  setPtr->insert(facei);
1654  }
1655 
1656  break;
1657  }
1658 
1659  prevN = triN;
1660  }
1661  else if (minTwist > 0)
1662  {
1663  nWarped++;
1664 
1665  if (setPtr)
1666  {
1667  setPtr->insert(facei);
1668  }
1669 
1670  break;
1671  }
1672  }
1673  while (fp != startFp);
1674  }
1675  }
1676  }
1677 
1678 
1679  reduce(nWarped, sumOp<label>());
1680 
1681  if (report)
1682  {
1683  if (nWarped> 0)
1684  {
1685  Info<< "There are " << nWarped
1686  << " faces with cosine of the angle"
1687  << " between consecutive triangle normals less than "
1688  << minTwist << nl << endl;
1689  }
1690  else
1691  {
1692  Info<< "All faces are flat in that the cosine of the angle"
1693  << " between consecutive triangle normals is less than "
1694  << minTwist << nl << endl;
1695  }
1696  }
1697 
1698  if (nWarped > 0)
1699  {
1700  if (report)
1701  {
1703  << nWarped << " faces with severe warpage "
1704  << "(cosine of the angle between consecutive triangle normals"
1705  << " < " << minTwist << ") found.\n"
1706  << endl;
1707  }
1708 
1709  return true;
1710  }
1711  else
1712  {
1713  return false;
1714  }
1715 }
1716 
1717 
1719 (
1720  const bool report,
1721  const scalar minFlatness,
1722  const polyMesh& mesh,
1723  const vectorField& faceAreas,
1724  const vectorField& faceCentres,
1725  const pointField& p,
1726  const labelList& checkFaces,
1727  labelHashSet* setPtr
1728 )
1729 {
1730  if (minFlatness < -small || minFlatness > 1+small)
1731  {
1733  << "minFlatness should be [0..1] but is now " << minFlatness
1734  << abort(FatalError);
1735  }
1736 
1737  const faceList& fcs = mesh.faces();
1738 
1739  label nWarped = 0;
1740 
1741  forAll(checkFaces, i)
1742  {
1743  const label facei = checkFaces[i];
1744 
1745  const face& f = fcs[facei];
1746 
1747  if (f.size() > 3)
1748  {
1749  const point& fc = faceCentres[facei];
1750 
1751  // Sum triangle areas
1752  scalar sumArea = 0.0;
1753 
1754  forAll(f, fp)
1755  {
1756  sumArea += triPointRef
1757  (
1758  p[f[fp]],
1759  p[f.nextLabel(fp)],
1760  fc
1761  ).mag();
1762  }
1763 
1764  if (sumArea/mag(faceAreas[facei]) < minFlatness)
1765  {
1766  nWarped++;
1767 
1768  if (setPtr)
1769  {
1770  setPtr->insert(facei);
1771  }
1772  }
1773  }
1774  }
1775 
1776  reduce(nWarped, sumOp<label>());
1777 
1778  if (report)
1779  {
1780  if (nWarped> 0)
1781  {
1782  Info<< "There are " << nWarped
1783  << " faces with area of individual triangles"
1784  << " compared to overall area less than "
1785  << minFlatness << nl << endl;
1786  }
1787  else
1788  {
1789  Info<< "All faces are flat in that the area of individual triangles"
1790  << " compared to overall area is less than "
1791  << minFlatness << nl << endl;
1792  }
1793  }
1794 
1795  if (nWarped > 0)
1796  {
1797  if (report)
1798  {
1800  << nWarped << " non-flat faces "
1801  << "(area of individual triangles"
1802  << " compared to overall area"
1803  << " < " << minFlatness << ") found.\n"
1804  << endl;
1805  }
1806 
1807  return true;
1808  }
1809  else
1810  {
1811  return false;
1812  }
1813 }
1814 
1815 
1817 (
1818  const bool report,
1819  const scalar minArea,
1820  const polyMesh& mesh,
1821  const vectorField& faceAreas,
1822  const labelList& checkFaces,
1823  labelHashSet* setPtr
1824 )
1825 {
1826  label nZeroArea = 0;
1827 
1828  forAll(checkFaces, i)
1829  {
1830  const label facei = checkFaces[i];
1831 
1832  if (mag(faceAreas[facei]) < minArea)
1833  {
1834  if (setPtr)
1835  {
1836  setPtr->insert(facei);
1837  }
1838  nZeroArea++;
1839  }
1840  }
1841 
1842 
1843  reduce(nZeroArea, sumOp<label>());
1844 
1845  if (report)
1846  {
1847  if (nZeroArea > 0)
1848  {
1849  Info<< "There are " << nZeroArea
1850  << " faces with area < " << minArea << '.' << nl << endl;
1851  }
1852  else
1853  {
1854  Info<< "All faces have area > " << minArea << '.' << nl << endl;
1855  }
1856  }
1857 
1858  if (nZeroArea > 0)
1859  {
1860  if (report)
1861  {
1863  << nZeroArea << " faces with area < " << minArea
1864  << " found.\n"
1865  << endl;
1866  }
1867 
1868  return true;
1869  }
1870  else
1871  {
1872  return false;
1873  }
1874 }
1875 
1876 
1878 (
1879  const bool report,
1880  const scalar warnDet,
1881  const polyMesh& mesh,
1882  const vectorField& faceAreas,
1883  const labelList& checkFaces,
1884  labelHashSet* setPtr
1885 )
1886 {
1887  const cellList& cells = mesh.cells();
1888 
1889  scalar minDet = great;
1890  scalar sumDet = 0.0;
1891  label nSumDet = 0;
1892  label nWarnDet = 0;
1893 
1894  const labelList affectedCells(getAffectedCells(mesh, checkFaces));
1895 
1896  forAll(affectedCells, i)
1897  {
1898  const cell& cFaces = cells[affectedCells[i]];
1899 
1900  tensor areaSum(Zero);
1901  scalar magAreaSum = 0;
1902 
1903  forAll(cFaces, cFacei)
1904  {
1905  const label facei = cFaces[cFacei];
1906 
1907  const scalar magArea = mag(faceAreas[facei]);
1908 
1909  magAreaSum += magArea;
1910  areaSum += faceAreas[facei]*(faceAreas[facei]/(magArea + vSmall));
1911  }
1912 
1913  // Scalad by 1/(cube root(1/3)) = 27
1914  const scalar scaledDet = 27*det(areaSum/(magAreaSum + vSmall));
1915 
1916  minDet = min(minDet, scaledDet);
1917  sumDet += scaledDet;
1918  nSumDet++;
1919 
1920  if (scaledDet < warnDet)
1921  {
1922  if (setPtr)
1923  {
1924  // Insert all faces of the cell.
1925  forAll(cFaces, cFacei)
1926  {
1927  const label facei = cFaces[cFacei];
1928  setPtr->insert(facei);
1929  }
1930  }
1931  nWarnDet++;
1932  }
1933  }
1934 
1935  reduce(minDet, minOp<scalar>());
1936  reduce(sumDet, sumOp<scalar>());
1937  reduce(nSumDet, sumOp<label>());
1938  reduce(nWarnDet, sumOp<label>());
1939 
1940  if (report)
1941  {
1942  if (nSumDet > 0)
1943  {
1944  Info<< "Cell determinant (1 = uniform cube) : average = "
1945  << sumDet / nSumDet << " min = " << minDet << endl;
1946  }
1947 
1948  if (nWarnDet > 0)
1949  {
1950  Info<< "There are " << nWarnDet
1951  << " cells with determinant < " << warnDet << '.' << nl
1952  << endl;
1953  }
1954  else
1955  {
1956  Info<< "All faces have determinant > " << warnDet << '.' << nl
1957  << endl;
1958  }
1959  }
1960 
1961  if (nWarnDet > 0)
1962  {
1963  if (report)
1964  {
1966  << nWarnDet << " cells with determinant < " << warnDet
1967  << " found.\n"
1968  << endl;
1969  }
1970 
1971  return true;
1972  }
1973  else
1974  {
1975  return false;
1976  }
1977 }
1978 
1979 
1980 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:227
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
T & first()
Return the first element of the list.
Definition: UListI.H:114
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1345
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
virtual const faceList & faces() const =0
Return faces.
virtual const labelList & faceOwner() const =0
Face face-owner addressing.
label nInternalFaces() const
label nFaces() const
virtual const labelList & faceNeighbour() const =0
Face face-neighbour addressing.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
const cellList & cells() const
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:80
vector area() const
Return vector area.
Definition: triangleI.H:96
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const pointField & points
const cellShapeList & cells
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
bool checkFacePyramids(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
Check face pyramid volumes.
tmp< scalarField > faceSkewness(const polyMesh &mesh, const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate skewness field.
Definition: polyMeshCheck.C:89
bool checkTriangleTwist(const bool report, const scalar minTwist, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Check consecutive face-triangle (from face-centre decomposition) normals.
bool checkFaceTwist(const bool report, const scalar minTwist, const polyMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Check the difference between normals of individual face-triangles (from.
labelList getAffectedCells(const primitiveMesh &mesh, const labelList &changedFaces)
bool checkFaceArea(const bool report, const scalar minArea, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Check for small faces.
bool checkFaceWeights(const bool report, const scalar warnWeight, const polyMesh &mesh, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
Check interpolation weights (0.5 for regular mesh)
bool checkCellDeterminant(const polyMesh &mesh, const bool report, labelHashSet *setPtr)
bool checkFaceAngles(const bool report, const scalar maxConcave, const polyMesh &mesh, const vectorField &faceAreas, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Check convexity of angles in a face. See primitiveMesh for explanation.
bool checkVolRatio(const polyMesh &mesh, const bool report, const scalar minRatio=0.01, labelHashSet *setPtr=nullptr)
Check for neighbouring cell volumes.
scalar boundaryFaceSkewness(const primitiveMesh &mesh, const pointField &p, const vectorField &fCtrs, const vectorField &fAreas, const label facei, const point &ownCc)
Skewness of single boundary face.
scalar checkNonOrtho(const primitiveMesh &mesh, const bool report, const scalar severeNonorthogonalityThreshold, const label facei, const vector &s, const vector &d, label &severeNonOrth, label &errorNonOrth, labelHashSet *setPtr)
bool checkFaceTets(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
Check face tetrahedra volumes.
bool checkFaceFlatness(const bool report, const scalar minFlatness, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Check for face areas v.s. sum of face-triangle (from face-centre.
bool checkFaceOrthogonality(const polyMesh &mesh, const scalar nonOrthThreshold, const bool report=false, labelHashSet *setPtr=nullptr)
Check face orthogonality.
bool checkFaceSkewness(const polyMesh &mesh, const scalar skewThreshold, const bool report=false, labelHashSet *setPtr=nullptr)
Check face skewness.
bool checkFaceTet(const primitiveMesh &mesh, const bool report, const scalar minTetQuality, const pointField &p, const label facei, const point &fc, const point &cc, labelHashSet *setPtr)
Namespace for OpenFOAM.
pyramid< point, const point &, const face & > pyramidPointFaceRef
dimensionedScalar det(const dimensionedSphericalTensor &dt)
static const zero Zero
Definition: zero.H:97
dimensionedScalar asin(const dimensionedScalar &ds)
scalar radToDeg(const scalar rad)
Convert radians to degrees.
scalar degToRad(const scalar deg)
Convert degrees to radians.
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:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar sin(const dimensionedScalar &ds)
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)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
static const char nl
Definition: Ostream.H:266
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
tetrahedron< point, const point & > tetPointRef
Definition: tetPointRef.H:44
labelList f(nPoints)
volScalarField & p