primitiveMeshCheck.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) 2011-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 "primitiveMesh.H"
27 #include "pyramidPointFaceRef.H"
28 #include "ListOps.H"
29 #include "unitConversion.H"
30 #include "SortableList.H"
31 #include "EdgeMap.H"
32 #include "primitiveMeshTools.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 Foam::scalar Foam::primitiveMesh::closedThreshold_ = 1.0e-6;
37 Foam::scalar Foam::primitiveMesh::aspectThreshold_ = 1000;
38 Foam::scalar Foam::primitiveMesh::nonOrthThreshold_ = 70; // deg
40 Foam::scalar Foam::primitiveMesh::planarCosAngle_ = 1.0e-6;
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 (
47  const vectorField& areas,
48  const bool report,
49  const PackedBoolList& internalOrCoupledFaces
50 ) const
51 {
52  if (debug)
53  {
54  Info<< "bool primitiveMesh::checkClosedBoundary("
55  << "const bool) const: "
56  << "checking whether the boundary is closed" << endl;
57  }
58 
59  // Loop through all boundary faces and sum up the face area vectors.
60  // For a closed boundary, this should be zero in all vector components
61 
62  vector sumClosed(vector::zero);
63  scalar sumMagClosedBoundary = 0;
64 
65  for (label faceI = nInternalFaces(); faceI < areas.size(); faceI++)
66  {
67  if (!internalOrCoupledFaces.size() || !internalOrCoupledFaces[faceI])
68  {
69  sumClosed += areas[faceI];
70  sumMagClosedBoundary += mag(areas[faceI]);
71  }
72  }
73 
74  reduce(sumClosed, sumOp<vector>());
75  reduce(sumMagClosedBoundary, sumOp<scalar>());
76 
77  vector openness = sumClosed/(sumMagClosedBoundary + VSMALL);
78 
79  if (cmptMax(cmptMag(openness)) > closedThreshold_)
80  {
81  if (debug || report)
82  {
83  Info<< " ***Boundary openness " << openness
84  << " possible hole in boundary description."
85  << endl;
86  }
87 
88  return true;
89  }
90  else
91  {
92  if (debug || report)
93  {
94  Info<< " Boundary openness " << openness << " OK."
95  << endl;
96  }
97 
98  return false;
99  }
100 }
101 
102 
104 (
105  const vectorField& faceAreas,
106  const scalarField& cellVolumes,
107  const bool report,
108  labelHashSet* setPtr,
109  labelHashSet* aspectSetPtr,
110  const Vector<label>& meshD
111 ) const
112 {
113  if (debug)
114  {
115  Info<< "bool primitiveMesh::checkClosedCells("
116  << "const bool, labelHashSet*, labelHashSet*"
117  << ", const Vector<label>&) const: "
118  << "checking whether cells are closed" << endl;
119  }
120 
121  // Check that all cells labels are valid
122  const cellList& c = cells();
123 
124  label nErrorClosed = 0;
125 
126  forAll(c, cI)
127  {
128  const cell& curCell = c[cI];
129 
130  if (min(curCell) < 0 || max(curCell) > nFaces())
131  {
132  if (setPtr)
133  {
134  setPtr->insert(cI);
135  }
136 
137  nErrorClosed++;
138  }
139  }
140 
141  if (nErrorClosed > 0)
142  {
143  if (debug || report)
144  {
145  Info<< " ***Cells with invalid face labels found, number of cells "
146  << nErrorClosed << endl;
147  }
148 
149  return true;
150  }
151 
152 
153  scalarField openness;
154  scalarField aspectRatio;
155  primitiveMeshTools::cellClosedness
156  (
157  *this,
158  meshD,
159  faceAreas,
160  cellVolumes,
161  openness,
162  aspectRatio
163  );
164 
165  label nOpen = 0;
166  scalar maxOpennessCell = max(openness);
167  label nAspect = 0;
168  scalar maxAspectRatio = max(aspectRatio);
169 
170  // Check the sums
171  forAll(openness, cellI)
172  {
173  if (openness[cellI] > closedThreshold_)
174  {
175  if (setPtr)
176  {
177  setPtr->insert(cellI);
178  }
179 
180  nOpen++;
181  }
182 
183  if (aspectRatio[cellI] > aspectThreshold_)
184  {
185  if (aspectSetPtr)
186  {
187  aspectSetPtr->insert(cellI);
188  }
189 
190  nAspect++;
191  }
192  }
193 
194  reduce(nOpen, sumOp<label>());
195  reduce(maxOpennessCell, maxOp<scalar>());
196 
197  reduce(nAspect, sumOp<label>());
198  reduce(maxAspectRatio, maxOp<scalar>());
199 
200 
201  if (nOpen > 0)
202  {
203  if (debug || report)
204  {
205  Info<< " ***Open cells found, max cell openness: "
206  << maxOpennessCell << ", number of open cells " << nOpen
207  << endl;
208  }
209 
210  return true;
211  }
212 
213  if (nAspect > 0)
214  {
215  if (debug || report)
216  {
217  Info<< " ***High aspect ratio cells found, Max aspect ratio: "
218  << maxAspectRatio
219  << ", number of cells " << nAspect
220  << endl;
221  }
222 
223  return true;
224  }
225 
226  if (debug || report)
227  {
228  Info<< " Max cell openness = " << maxOpennessCell << " OK." << nl
229  << " Max aspect ratio = " << maxAspectRatio << " OK."
230  << endl;
231  }
232 
233  return false;
234 }
235 
236 
238 (
239  const vectorField& faceAreas,
240  const bool report,
241  const bool detailedReport,
242  labelHashSet* setPtr
243 ) const
244 {
245  if (debug)
246  {
247  Info<< "bool primitiveMesh::checkFaceAreas("
248  << "const bool, labelHashSet*) const: "
249  << "checking face area magnitudes" << endl;
250  }
251 
252  const scalarField magFaceAreas(mag(faceAreas));
253 
254  scalar minArea = GREAT;
255  scalar maxArea = -GREAT;
256 
257  forAll(magFaceAreas, faceI)
258  {
259  if (magFaceAreas[faceI] < VSMALL)
260  {
261  if (setPtr)
262  {
263  setPtr->insert(faceI);
264  }
265  if (detailedReport)
266  {
267  if (isInternalFace(faceI))
268  {
269  Pout<< "Zero or negative face area detected for "
270  << "internal face "<< faceI << " between cells "
271  << faceOwner()[faceI] << " and "
272  << faceNeighbour()[faceI]
273  << ". Face area magnitude = " << magFaceAreas[faceI]
274  << endl;
275  }
276  else
277  {
278  Pout<< "Zero or negative face area detected for "
279  << "boundary face " << faceI << " next to cell "
280  << faceOwner()[faceI] << ". Face area magnitude = "
281  << magFaceAreas[faceI] << endl;
282  }
283  }
284  }
285 
286  minArea = min(minArea, magFaceAreas[faceI]);
287  maxArea = max(maxArea, magFaceAreas[faceI]);
288  }
289 
290  reduce(minArea, minOp<scalar>());
291  reduce(maxArea, maxOp<scalar>());
292 
293  if (minArea < VSMALL)
294  {
295  if (debug || report)
296  {
297  Info<< " ***Zero or negative face area detected. "
298  "Minimum area: " << minArea << endl;
299  }
300 
301  return true;
302  }
303  else
304  {
305  if (debug || report)
306  {
307  Info<< " Minimum face area = " << minArea
308  << ". Maximum face area = " << maxArea
309  << ". Face area magnitudes OK." << endl;
310  }
311 
312  return false;
313  }
314 }
315 
316 
318 (
319  const scalarField& vols,
320  const bool report,
321  const bool detailedReport,
322  labelHashSet* setPtr
323 ) const
324 {
325  if (debug)
326  {
327  Info<< "bool primitiveMesh::checkCellVolumes("
328  << "const bool, labelHashSet*) const: "
329  << "checking cell volumes" << endl;
330  }
331 
332  scalar minVolume = GREAT;
333  scalar maxVolume = -GREAT;
334 
335  label nNegVolCells = 0;
336 
337  forAll(vols, cellI)
338  {
339  if (vols[cellI] < VSMALL)
340  {
341  if (setPtr)
342  {
343  setPtr->insert(cellI);
344  }
345  if (detailedReport)
346  {
347  Pout<< "Zero or negative cell volume detected for cell "
348  << cellI << ". Volume = " << vols[cellI] << endl;
349  }
350 
351  nNegVolCells++;
352  }
353 
354  minVolume = min(minVolume, vols[cellI]);
355  maxVolume = max(maxVolume, vols[cellI]);
356  }
357 
358  reduce(minVolume, minOp<scalar>());
359  reduce(maxVolume, maxOp<scalar>());
360  reduce(nNegVolCells, sumOp<label>());
361 
362  if (minVolume < VSMALL)
363  {
364  if (debug || report)
365  {
366  Info<< " ***Zero or negative cell volume detected. "
367  << "Minimum negative volume: " << minVolume
368  << ", Number of negative volume cells: " << nNegVolCells
369  << endl;
370  }
371 
372  return true;
373  }
374  else
375  {
376  if (debug || report)
377  {
378  Info<< " Min volume = " << minVolume
379  << ". Max volume = " << maxVolume
380  << ". Total volume = " << gSum(vols)
381  << ". Cell volumes OK." << endl;
382  }
383 
384  return false;
385  }
386 }
387 
388 
390 (
391  const vectorField& fAreas,
392  const vectorField& cellCtrs,
393  const bool report,
394  labelHashSet* setPtr
395 ) const
396 {
397  if (debug)
398  {
399  Info<< "bool primitiveMesh::checkFaceOrthogonality("
400  << "const bool, labelHashSet*) const: "
401  << "checking mesh non-orthogonality" << endl;
402  }
403 
404 
405  tmp<scalarField> tortho = primitiveMeshTools::faceOrthogonality
406  (
407  *this,
408  fAreas,
409  cellCtrs
410  );
411  const scalarField& ortho = tortho();
412 
413  // Severe nonorthogonality threshold
414  const scalar severeNonorthogonalityThreshold =
415  ::cos(degToRad(nonOrthThreshold_));
416 
417  scalar minDDotS = min(ortho);
418 
419  scalar sumDDotS = sum(ortho);
420 
421  label severeNonOrth = 0;
422 
423  label errorNonOrth = 0;
424 
425 
426  forAll(ortho, faceI)
427  {
428  if (ortho[faceI] < severeNonorthogonalityThreshold)
429  {
430  if (ortho[faceI] > SMALL)
431  {
432  if (setPtr)
433  {
434  setPtr->insert(faceI);
435  }
436 
437  severeNonOrth++;
438  }
439  else
440  {
441  if (setPtr)
442  {
443  setPtr->insert(faceI);
444  }
445 
446  errorNonOrth++;
447  }
448  }
449  }
450 
451  reduce(minDDotS, minOp<scalar>());
452  reduce(sumDDotS, sumOp<scalar>());
453  reduce(severeNonOrth, sumOp<label>());
454  reduce(errorNonOrth, sumOp<label>());
455 
456  if (debug || report)
457  {
458  label neiSize = ortho.size();
459  reduce(neiSize, sumOp<label>());
460 
461  if (neiSize > 0)
462  {
463  if (debug || report)
464  {
465  Info<< " Mesh non-orthogonality Max: "
466  << radToDeg(::acos(minDDotS))
467  << " average: " << radToDeg(::acos(sumDDotS/neiSize))
468  << endl;
469  }
470  }
471 
472  if (severeNonOrth > 0)
473  {
474  Info<< " *Number of severely non-orthogonal faces: "
475  << severeNonOrth << "." << endl;
476  }
477  }
478 
479  if (errorNonOrth > 0)
480  {
481  if (debug || report)
482  {
483  Info<< " ***Number of non-orthogonality errors: "
484  << errorNonOrth << "." << endl;
485  }
486 
487  return true;
488  }
489  else
490  {
491  if (debug || report)
492  {
493  Info<< " Non-orthogonality check OK." << endl;
494  }
495 
496  return false;
497  }
498 }
499 
500 
502 (
503  const pointField& points,
504  const vectorField& ctrs,
505  const bool report,
506  const bool detailedReport,
507  const scalar minPyrVol,
508  labelHashSet* setPtr
509 ) const
510 {
511  if (debug)
512  {
513  Info<< "bool primitiveMesh::checkFacePyramids("
514  << "const bool, const scalar, labelHashSet*) const: "
515  << "checking face orientation" << endl;
516  }
517 
518  const labelList& own = faceOwner();
519  const labelList& nei = faceNeighbour();
520  const faceList& f = faces();
521 
522 
523  scalarField ownPyrVol;
524  scalarField neiPyrVol;
525  primitiveMeshTools::facePyramidVolume
526  (
527  *this,
528  points,
529  ctrs,
530  ownPyrVol,
531  neiPyrVol
532  );
533 
534 
535  label nErrorPyrs = 0;
536 
537  forAll(ownPyrVol, faceI)
538  {
539  if (ownPyrVol[faceI] < minPyrVol)
540  {
541  if (setPtr)
542  {
543  setPtr->insert(faceI);
544  }
545  if (detailedReport)
546  {
547  Pout<< "Negative pyramid volume: " << ownPyrVol[faceI]
548  << " for face " << faceI << " " << f[faceI]
549  << " and owner cell: " << own[faceI] << endl
550  << "Owner cell vertex labels: "
551  << cells()[own[faceI]].labels(faces())
552  << endl;
553  }
554 
555  nErrorPyrs++;
556  }
557 
558  if (isInternalFace(faceI))
559  {
560  if (neiPyrVol[faceI] < minPyrVol)
561  {
562  if (setPtr)
563  {
564  setPtr->insert(faceI);
565  }
566  if (detailedReport)
567  {
568  Pout<< "Negative pyramid volume: " << neiPyrVol[faceI]
569  << " for face " << faceI << " " << f[faceI]
570  << " and neighbour cell: " << nei[faceI] << nl
571  << "Neighbour cell vertex labels: "
572  << cells()[nei[faceI]].labels(faces())
573  << endl;
574  }
575  nErrorPyrs++;
576  }
577  }
578  }
579 
580  reduce(nErrorPyrs, sumOp<label>());
581 
582  if (nErrorPyrs > 0)
583  {
584  if (debug || report)
585  {
586  Info<< " ***Error in face pyramids: "
587  << nErrorPyrs << " faces are incorrectly oriented."
588  << endl;
589  }
590 
591  return true;
592  }
593  else
594  {
595  if (debug || report)
596  {
597  Info<< " Face pyramids OK." << endl;
598  }
599 
600  return false;
601  }
602 }
603 
604 
606 (
607  const pointField& points,
608  const vectorField& fCtrs,
609  const vectorField& fAreas,
610  const vectorField& cellCtrs,
611  const bool report,
612  labelHashSet* setPtr
613 ) const
614 {
615  if (debug)
616  {
617  Info<< "bool primitiveMesh::checkFaceSkewnesss("
618  << "const bool, labelHashSet*) const: "
619  << "checking face skewness" << endl;
620  }
621 
622  // Warn if the skew correction vector is more than skewWarning times
623  // larger than the face area vector
624 
625  tmp<scalarField> tskewness = primitiveMeshTools::faceSkewness
626  (
627  *this,
628  points,
629  fCtrs,
630  fAreas,
631  cellCtrs
632  );
633  const scalarField& skewness = tskewness();
634 
635  scalar maxSkew = max(skewness);
636  label nWarnSkew = 0;
637 
638  forAll(skewness, faceI)
639  {
640  // Check if the skewness vector is greater than the PN vector.
641  // This does not cause trouble but is a good indication of a poor mesh.
642  if (skewness[faceI] > skewThreshold_)
643  {
644  if (setPtr)
645  {
646  setPtr->insert(faceI);
647  }
648 
649  nWarnSkew++;
650  }
651  }
652 
653  reduce(maxSkew, maxOp<scalar>());
654  reduce(nWarnSkew, sumOp<label>());
655 
656  if (nWarnSkew > 0)
657  {
658  if (debug || report)
659  {
660  Info<< " ***Max skewness = " << maxSkew
661  << ", " << nWarnSkew << " highly skew faces detected"
662  " which may impair the quality of the results"
663  << endl;
664  }
665 
666  return true;
667  }
668  else
669  {
670  if (debug || report)
671  {
672  Info<< " Max skewness = " << maxSkew << " OK." << endl;
673  }
674 
675  return false;
676  }
677 }
678 
679 
680 // Check convexity of angles in a face. Allow a slight non-convexity.
681 // E.g. maxDeg = 10 allows for angles < 190 (or 10 degrees concavity)
682 // (if truly concave and points not visible from face centre the face-pyramid
683 // check in checkMesh will fail)
685 (
686  const pointField& points,
687  const vectorField& faceAreas,
688  const bool report,
689  const scalar maxDeg,
690  labelHashSet* setPtr
691 ) const
692 {
693  if (debug)
694  {
695  Info<< "bool primitiveMesh::checkFaceAngles"
696  << "(const bool, const scalar, labelHashSet*) const: "
697  << "checking face angles" << endl;
698  }
699 
700  if (maxDeg < -SMALL || maxDeg > 180+SMALL)
701  {
703  << "maxDeg should be [0..180] but is now " << maxDeg
704  << exit(FatalError);
705  }
706 
707  const scalar maxSin = Foam::sin(degToRad(maxDeg));
708 
709 
710  tmp<scalarField> tfaceAngles = primitiveMeshTools::faceConcavity
711  (
712  maxSin,
713  *this,
714  points,
715  faceAreas
716  );
717  const scalarField& faceAngles = tfaceAngles();
718 
719  scalar maxEdgeSin = max(faceAngles);
720 
721  label nConcave = 0;
722 
723  forAll(faceAngles, faceI)
724  {
725  if (faceAngles[faceI] > SMALL)
726  {
727  nConcave++;
728 
729  if (setPtr)
730  {
731  setPtr->insert(faceI);
732  }
733  }
734  }
735 
736  reduce(nConcave, sumOp<label>());
737  reduce(maxEdgeSin, maxOp<scalar>());
738 
739  if (nConcave > 0)
740  {
741  scalar maxConcaveDegr =
742  radToDeg(Foam::asin(Foam::min(1.0, maxEdgeSin)));
743 
744  if (debug || report)
745  {
746  Info<< " *There are " << nConcave
747  << " faces with concave angles between consecutive"
748  << " edges. Max concave angle = " << maxConcaveDegr
749  << " degrees." << endl;
750  }
751 
752  return true;
753  }
754  else
755  {
756  if (debug || report)
757  {
758  Info<< " All angles in faces OK." << endl;
759  }
760 
761  return false;
762  }
763 }
764 
765 
767 (
768  const pointField& points,
769  const vectorField& faceCentres,
770  const vectorField& faceAreas,
771  const bool report,
772  const scalar warnFlatness,
773  labelHashSet* setPtr
774 ) const
775 {
776  if (debug)
777  {
778  Info<< "bool primitiveMesh::checkFaceFlatness"
779  << "(const bool, const scalar, labelHashSet*) const: "
780  << "checking face flatness" << endl;
781  }
782 
783  if (warnFlatness < 0 || warnFlatness > 1)
784  {
786  << "warnFlatness should be [0..1] but is now " << warnFlatness
787  << exit(FatalError);
788  }
789 
790  const faceList& fcs = faces();
791 
792  tmp<scalarField> tfaceFlatness = primitiveMeshTools::faceFlatness
793  (
794  *this,
795  points,
796  faceCentres,
797  faceAreas
798  );
799  const scalarField& faceFlatness = tfaceFlatness();
800 
801  scalarField magAreas(mag(faceAreas));
802 
803  scalar minFlatness = GREAT;
804  scalar sumFlatness = 0;
805  label nSummed = 0;
806  label nWarped = 0;
807 
808  forAll(faceFlatness, faceI)
809  {
810  if (fcs[faceI].size() > 3 && magAreas[faceI] > VSMALL)
811  {
812  sumFlatness += faceFlatness[faceI];
813  nSummed++;
814 
815  minFlatness = min(minFlatness, faceFlatness[faceI]);
816 
817  if (faceFlatness[faceI] < warnFlatness)
818  {
819  nWarped++;
820 
821  if (setPtr)
822  {
823  setPtr->insert(faceI);
824  }
825  }
826  }
827  }
828 
829 
830  reduce(nWarped, sumOp<label>());
831  reduce(minFlatness, minOp<scalar>());
832 
833  reduce(nSummed, sumOp<label>());
834  reduce(sumFlatness, sumOp<scalar>());
835 
836  if (debug || report)
837  {
838  if (nSummed > 0)
839  {
840  Info<< " Face flatness (1 = flat, 0 = butterfly) : min = "
841  << minFlatness << " average = " << sumFlatness / nSummed
842  << endl;
843  }
844  }
845 
846 
847  if (nWarped> 0)
848  {
849  if (debug || report)
850  {
851  Info<< " *There are " << nWarped
852  << " faces with ratio between projected and actual area < "
853  << warnFlatness << endl;
854 
855  Info<< " Minimum ratio (minimum flatness, maximum warpage) = "
856  << minFlatness << endl;
857  }
858 
859  return true;
860  }
861  else
862  {
863  if (debug || report)
864  {
865  Info<< " All face flatness OK." << endl;
866  }
867 
868  return false;
869  }
870 }
871 
872 
874 (
875  const vectorField& fAreas,
876  const pointField& fCentres,
877  const bool report,
878  labelHashSet* setPtr
879 ) const
880 {
881  if (debug)
882  {
883  Info<< "bool primitiveMesh::checkConcaveCells(const bool"
884  << ", labelHashSet*) const: "
885  << "checking for concave cells" << endl;
886  }
887 
888  const cellList& c = cells();
889  const labelList& fOwner = faceOwner();
890 
891  label nConcaveCells = 0;
892 
893  forAll(c, cellI)
894  {
895  const cell& cFaces = c[cellI];
896 
897  bool concave = false;
898 
899  forAll(cFaces, i)
900  {
901  if (concave)
902  {
903  break;
904  }
905 
906  label fI = cFaces[i];
907 
908  const point& fC = fCentres[fI];
909 
910  vector fN = fAreas[fI];
911 
912  fN /= max(mag(fN), VSMALL);
913 
914  // Flip normal if required so that it is always pointing out of
915  // the cell
916  if (fOwner[fI] != cellI)
917  {
918  fN *= -1;
919  }
920 
921  // Is the centre of any other face of the cell on the
922  // wrong side of the plane of this face?
923 
924  forAll(cFaces, j)
925  {
926  if (j != i)
927  {
928  label fJ = cFaces[j];
929 
930  const point& pt = fCentres[fJ];
931 
932  // If the cell is concave, the point will be on the
933  // positive normal side of the plane of f, defined by
934  // its centre and normal, and the angle between (pt -
935  // fC) and fN will be less than 90 degrees, so the dot
936  // product will be positive.
937 
938  vector pC = (pt - fC);
939 
940  pC /= max(mag(pC), VSMALL);
941 
942  if ((pC & fN) > -planarCosAngle_)
943  {
944  // Concave or planar face
945 
946  concave = true;
947 
948  if (setPtr)
949  {
950  setPtr->insert(cellI);
951  }
952 
953  nConcaveCells++;
954 
955  break;
956  }
957  }
958  }
959  }
960  }
961 
962  reduce(nConcaveCells, sumOp<label>());
963 
964  if (nConcaveCells > 0)
965  {
966  if (debug || report)
967  {
968  Info<< " ***Concave cells (using face planes) found,"
969  << " number of cells: " << nConcaveCells << endl;
970  }
971 
972  return true;
973  }
974  else
975  {
976  if (debug || report)
977  {
978  Info<< " Concave cell check OK." << endl;
979  }
980 
981  return false;
982  }
983 
984  return false;
985 }
986 
987 
988 // Topological tests
989 
991 (
992  const bool report,
993  labelHashSet* setPtr
994 ) const
995 {
996  if (debug)
997  {
998  Info<< "bool primitiveMesh::checkUpperTriangular("
999  << "const bool, labelHashSet*) const: "
1000  << "checking face ordering" << endl;
1001  }
1002 
1003  // Check whether internal faces are ordered in the upper triangular order
1004  const labelList& own = faceOwner();
1005  const labelList& nei = faceNeighbour();
1006 
1007  const cellList& c = cells();
1008 
1009  label internal = nInternalFaces();
1010 
1011  // Has error occurred?
1012  bool error = false;
1013  // Have multiple faces been detected?
1014  label nMultipleCells = false;
1015 
1016  // Loop through faceCells once more and make sure that for internal cell
1017  // the first label is smaller
1018  for (label faceI = 0; faceI < internal; faceI++)
1019  {
1020  if (own[faceI] >= nei[faceI])
1021  {
1022  error = true;
1023 
1024  if (setPtr)
1025  {
1026  setPtr->insert(faceI);
1027  }
1028  }
1029  }
1030 
1031  // Loop through all cells. For each cell, find the face that is internal
1032  // and add it to the check list (upper triangular order).
1033  // Once the list is completed, check it against the faceCell list
1034 
1035  forAll(c, cellI)
1036  {
1037  const labelList& curFaces = c[cellI];
1038 
1039  // Neighbouring cells
1040  SortableList<label> nbr(curFaces.size());
1041 
1042  forAll(curFaces, i)
1043  {
1044  label faceI = curFaces[i];
1045 
1046  if (faceI >= nInternalFaces())
1047  {
1048  // Sort last
1049  nbr[i] = labelMax;
1050  }
1051  else
1052  {
1053  label nbrCellI = nei[faceI];
1054 
1055  if (nbrCellI == cellI)
1056  {
1057  nbrCellI = own[faceI];
1058  }
1059 
1060  if (cellI < nbrCellI)
1061  {
1062  // cellI is master
1063  nbr[i] = nbrCellI;
1064  }
1065  else
1066  {
1067  // nbrCell is master. Let it handle this face.
1068  nbr[i] = labelMax;
1069  }
1070  }
1071  }
1072 
1073  nbr.sort();
1074 
1075  // Now nbr holds the cellCells in incremental order. Check:
1076  // - neighbouring cells appear only once. Since nbr is sorted this
1077  // is simple check on consecutive elements
1078  // - faces indexed in same order as nbr are incrementing as well.
1079 
1080  label prevCell = nbr[0];
1081  label prevFace = curFaces[nbr.indices()[0]];
1082 
1083  bool hasMultipleFaces = false;
1084 
1085  for (label i = 1; i < nbr.size(); i++)
1086  {
1087  label thisCell = nbr[i];
1088  label thisFace = curFaces[nbr.indices()[i]];
1089 
1090  if (thisCell == labelMax)
1091  {
1092  break;
1093  }
1094 
1095  if (thisCell == prevCell)
1096  {
1097  hasMultipleFaces = true;
1098 
1099  if (setPtr)
1100  {
1101  setPtr->insert(prevFace);
1102  setPtr->insert(thisFace);
1103  }
1104  }
1105  else if (thisFace < prevFace)
1106  {
1107  error = true;
1108 
1109  if (setPtr)
1110  {
1111  setPtr->insert(thisFace);
1112  }
1113  }
1114 
1115  prevCell = thisCell;
1116  prevFace = thisFace;
1117  }
1118 
1119  if (hasMultipleFaces)
1120  {
1121  nMultipleCells++;
1122  }
1123  }
1124 
1125  reduce(error, orOp<bool>());
1126  reduce(nMultipleCells, sumOp<label>());
1127 
1128  if ((debug || report) && nMultipleCells > 0)
1129  {
1130  Info<< " <<Found " << nMultipleCells
1131  << " neighbouring cells with multiple inbetween faces." << endl;
1132  }
1133 
1134  if (error)
1135  {
1136  if (debug || report)
1137  {
1138  Info<< " ***Faces not in upper triangular order." << endl;
1139  }
1140 
1141  return true;
1142  }
1143  else
1144  {
1145  if (debug || report)
1146  {
1147  Info<< " Upper triangular ordering OK." << endl;
1148  }
1149 
1150  return false;
1151  }
1152 }
1153 
1154 
1156 (
1157  const bool report,
1158  labelHashSet* setPtr
1159 ) const
1160 {
1161  if (debug)
1162  {
1163  Info<< "bool primitiveMesh::checkCellsZipUp("
1164  << "const bool, labelHashSet*) const: "
1165  << "checking topological cell openness" << endl;
1166  }
1167 
1168  label nOpenCells = 0;
1169 
1170  const faceList& f = faces();
1171  const cellList& c = cells();
1172 
1173  forAll(c, cellI)
1174  {
1175  const labelList& curFaces = c[cellI];
1176 
1177  const edgeList cellEdges = c[cellI].edges(f);
1178 
1179  labelList edgeUsage(cellEdges.size(), 0);
1180 
1181  forAll(curFaces, faceI)
1182  {
1183  edgeList curFaceEdges = f[curFaces[faceI]].edges();
1184 
1185  forAll(curFaceEdges, faceEdgeI)
1186  {
1187  const edge& curEdge = curFaceEdges[faceEdgeI];
1188 
1189  forAll(cellEdges, cellEdgeI)
1190  {
1191  if (cellEdges[cellEdgeI] == curEdge)
1192  {
1193  edgeUsage[cellEdgeI]++;
1194  break;
1195  }
1196  }
1197  }
1198  }
1199 
1200  edgeList singleEdges(cellEdges.size());
1201  label nSingleEdges = 0;
1202 
1203  forAll(edgeUsage, edgeI)
1204  {
1205  if (edgeUsage[edgeI] == 1)
1206  {
1207  singleEdges[nSingleEdges] = cellEdges[edgeI];
1208  nSingleEdges++;
1209  }
1210  else if (edgeUsage[edgeI] != 2)
1211  {
1212  if (setPtr)
1213  {
1214  setPtr->insert(cellI);
1215  }
1216  }
1217  }
1218 
1219  if (nSingleEdges > 0)
1220  {
1221  if (setPtr)
1222  {
1223  setPtr->insert(cellI);
1224  }
1225 
1226  nOpenCells++;
1227  }
1228  }
1229 
1230  reduce(nOpenCells, sumOp<label>());
1231 
1232  if (nOpenCells > 0)
1233  {
1234  if (debug || report)
1235  {
1236  Info<< " ***Open cells found, number of cells: " << nOpenCells
1237  << ". This problem may be fixable using the zipUpMesh utility."
1238  << endl;
1239  }
1240 
1241  return true;
1242  }
1243  else
1244  {
1245  if (debug || report)
1246  {
1247  Info<< " Topological cell zip-up check OK." << endl;
1248  }
1249 
1250  return false;
1251  }
1252 }
1253 
1254 
1255 // Vertices of face within point range and unique.
1257 (
1258  const bool report,
1259  labelHashSet* setPtr
1260 ) const
1261 {
1262  if (debug)
1263  {
1264  Info<< "bool primitiveMesh::checkFaceVertices("
1265  << "const bool, labelHashSet*) const: "
1266  << "checking face vertices" << endl;
1267  }
1268 
1269  // Check that all vertex labels are valid
1270  const faceList& f = faces();
1271 
1272  label nErrorFaces = 0;
1273 
1274  forAll(f, fI)
1275  {
1276  const face& curFace = f[fI];
1277 
1278  if (min(curFace) < 0 || max(curFace) > nPoints())
1279  {
1280  if (setPtr)
1281  {
1282  setPtr->insert(fI);
1283  }
1284 
1285  nErrorFaces++;
1286  }
1287 
1288  // Uniqueness of vertices
1289  labelHashSet facePoints(2*curFace.size());
1290 
1291  forAll(curFace, fp)
1292  {
1293  bool inserted = facePoints.insert(curFace[fp]);
1294 
1295  if (!inserted)
1296  {
1297  if (setPtr)
1298  {
1299  setPtr->insert(fI);
1300  }
1301 
1302  nErrorFaces++;
1303  }
1304  }
1305  }
1306 
1307  reduce(nErrorFaces, sumOp<label>());
1308 
1309  if (nErrorFaces > 0)
1310  {
1311  if (debug || report)
1312  {
1313  Info<< " Faces with invalid vertex labels found, "
1314  << " number of faces: " << nErrorFaces << endl;
1315  }
1316 
1317  return true;
1318  }
1319  else
1320  {
1321  if (debug || report)
1322  {
1323  Info<< " Face vertices OK." << endl;
1324  }
1325 
1326  return false;
1327  }
1328 }
1329 
1330 
1332 (
1333  const bool report,
1334  labelHashSet* setPtr
1335 ) const
1336 {
1337  if (debug)
1338  {
1339  Info<< "bool primitiveMesh::checkPoints"
1340  << "(const bool, labelHashSet*) const: "
1341  << "checking points" << endl;
1342  }
1343 
1344  label nFaceErrors = 0;
1345  label nCellErrors = 0;
1346 
1347  const labelListList& pf = pointFaces();
1348 
1349  forAll(pf, pointI)
1350  {
1351  if (pf[pointI].empty())
1352  {
1353  if (setPtr)
1354  {
1355  setPtr->insert(pointI);
1356  }
1357 
1358  nFaceErrors++;
1359  }
1360  }
1361 
1362 
1363  forAll(pf, pointI)
1364  {
1365  const labelList& pc = pointCells(pointI);
1366 
1367  if (pc.empty())
1368  {
1369  if (setPtr)
1370  {
1371  setPtr->insert(pointI);
1372  }
1373 
1374  nCellErrors++;
1375  }
1376  }
1377 
1378  reduce(nFaceErrors, sumOp<label>());
1379  reduce(nCellErrors, sumOp<label>());
1380 
1381  if (nFaceErrors > 0 || nCellErrors > 0)
1382  {
1383  if (debug || report)
1384  {
1385  Info<< " ***Unused points found in the mesh, "
1386  "number unused by faces: " << nFaceErrors
1387  << " number unused by cells: " << nCellErrors
1388  << endl;
1389  }
1390 
1391  return true;
1392  }
1393  else
1394  {
1395  if (debug || report)
1396  {
1397  Info<< " Point usage OK." << endl;
1398  }
1399 
1400  return false;
1401  }
1402 }
1403 
1404 
1405 // Check if all points on face are shared between faces.
1407 (
1408  const label faceI,
1409  const Map<label>& nCommonPoints,
1410  label& nBaffleFaces,
1411  labelHashSet* setPtr
1412 ) const
1413 {
1414  bool error = false;
1415 
1416  forAllConstIter(Map<label>, nCommonPoints, iter)
1417  {
1418  label nbFaceI = iter.key();
1419  label nCommon = iter();
1420 
1421  const face& curFace = faces()[faceI];
1422  const face& nbFace = faces()[nbFaceI];
1423 
1424  if (nCommon == nbFace.size() || nCommon == curFace.size())
1425  {
1426  if (nbFace.size() != curFace.size())
1427  {
1428  error = true;
1429  }
1430  else
1431  {
1432  nBaffleFaces++;
1433  }
1434 
1435  if (setPtr)
1436  {
1437  setPtr->insert(faceI);
1438  setPtr->insert(nbFaceI);
1439  }
1440  }
1441  }
1442 
1443  return error;
1444 }
1445 
1446 
1447 // Check that shared points are in consecutive order.
1449 (
1450  const label faceI,
1451  const Map<label>& nCommonPoints,
1452  labelHashSet* setPtr
1453 ) const
1454 {
1455  bool error = false;
1456 
1457  forAllConstIter(Map<label>, nCommonPoints, iter)
1458  {
1459  label nbFaceI = iter.key();
1460  label nCommon = iter();
1461 
1462  const face& curFace = faces()[faceI];
1463  const face& nbFace = faces()[nbFaceI];
1464 
1465  if
1466  (
1467  nCommon >= 2
1468  && nCommon != nbFace.size()
1469  && nCommon != curFace.size()
1470  )
1471  {
1472  forAll(curFace, fp)
1473  {
1474  // Get the index in the neighbouring face shared with curFace
1475  label nb = findIndex(nbFace, curFace[fp]);
1476 
1477  if (nb != -1)
1478  {
1479 
1480  // Check the whole face from nb onwards for shared vertices
1481  // with neighbouring face. Rule is that any shared vertices
1482  // should be consecutive on both faces i.e. if they are
1483  // vertices fp,fp+1,fp+2 on one face they should be
1484  // vertices nb, nb+1, nb+2 (or nb+2, nb+1, nb) on the
1485  // other face.
1486 
1487 
1488  // Vertices before and after on curFace
1489  label fpPlus1 = curFace.fcIndex(fp);
1490  label fpMin1 = curFace.rcIndex(fp);
1491 
1492  // Vertices before and after on nbFace
1493  label nbPlus1 = nbFace.fcIndex(nb);
1494  label nbMin1 = nbFace.rcIndex(nb);
1495 
1496  // Find order of walking by comparing next points on both
1497  // faces.
1498  label curInc = labelMax;
1499  label nbInc = labelMax;
1500 
1501  if (nbFace[nbPlus1] == curFace[fpPlus1])
1502  {
1503  curInc = 1;
1504  nbInc = 1;
1505  }
1506  else if (nbFace[nbPlus1] == curFace[fpMin1])
1507  {
1508  curInc = -1;
1509  nbInc = 1;
1510  }
1511  else if (nbFace[nbMin1] == curFace[fpMin1])
1512  {
1513  curInc = -1;
1514  nbInc = -1;
1515  }
1516  else
1517  {
1518  curInc = 1;
1519  nbInc = -1;
1520  }
1521 
1522 
1523  // Pass1: loop until start of common vertices found.
1524  label curNb = nb;
1525  label curFp = fp;
1526 
1527  do
1528  {
1529  curFp += curInc;
1530 
1531  if (curFp >= curFace.size())
1532  {
1533  curFp = 0;
1534  }
1535  else if (curFp < 0)
1536  {
1537  curFp = curFace.size()-1;
1538  }
1539 
1540  curNb += nbInc;
1541 
1542  if (curNb >= nbFace.size())
1543  {
1544  curNb = 0;
1545  }
1546  else if (curNb < 0)
1547  {
1548  curNb = nbFace.size()-1;
1549  }
1550  } while (curFace[curFp] == nbFace[curNb]);
1551 
1552 
1553  // Pass2: check equality walking from curFp, curNb
1554  // in opposite order.
1555 
1556  curInc = -curInc;
1557  nbInc = -nbInc;
1558 
1559  for (label commonI = 0; commonI < nCommon; commonI++)
1560  {
1561  curFp += curInc;
1562 
1563  if (curFp >= curFace.size())
1564  {
1565  curFp = 0;
1566  }
1567  else if (curFp < 0)
1568  {
1569  curFp = curFace.size()-1;
1570  }
1571 
1572  curNb += nbInc;
1573 
1574  if (curNb >= nbFace.size())
1575  {
1576  curNb = 0;
1577  }
1578  else if (curNb < 0)
1579  {
1580  curNb = nbFace.size()-1;
1581  }
1582 
1583  if (curFace[curFp] != nbFace[curNb])
1584  {
1585  if (setPtr)
1586  {
1587  setPtr->insert(faceI);
1588  setPtr->insert(nbFaceI);
1589  }
1590 
1591  error = true;
1592 
1593  break;
1594  }
1595  }
1596 
1597 
1598  // Done the curFace - nbFace combination.
1599  break;
1600  }
1601  }
1602  }
1603  }
1604 
1605  return error;
1606 }
1607 
1608 
1609 // Checks common vertices between faces. If more than 2 they should be
1610 // consecutive on both faces.
1612 (
1613  const bool report,
1614  labelHashSet* setPtr
1615 ) const
1616 {
1617  if (debug)
1618  {
1619  Info<< "bool primitiveMesh::checkFaceFaces(const bool, labelHashSet*)"
1620  << " const: " << "checking face-face connectivity" << endl;
1621  }
1622 
1623  const labelListList& pf = pointFaces();
1624 
1625  label nBaffleFaces = 0;
1626  label nErrorDuplicate = 0;
1627  label nErrorOrder = 0;
1628  Map<label> nCommonPoints(100);
1629 
1630  for (label faceI = 0; faceI < nFaces(); faceI++)
1631  {
1632  const face& curFace = faces()[faceI];
1633 
1634  // Calculate number of common points between current faceI and
1635  // neighbouring face. Store on map.
1636  nCommonPoints.clear();
1637 
1638  forAll(curFace, fp)
1639  {
1640  label pointI = curFace[fp];
1641 
1642  const labelList& nbs = pf[pointI];
1643 
1644  forAll(nbs, nbI)
1645  {
1646  label nbFaceI = nbs[nbI];
1647 
1648  if (faceI < nbFaceI)
1649  {
1650  // Only check once for each combination of two faces.
1651 
1652  Map<label>::iterator fnd = nCommonPoints.find(nbFaceI);
1653 
1654  if (fnd == nCommonPoints.end())
1655  {
1656  // First common vertex found.
1657  nCommonPoints.insert(nbFaceI, 1);
1658  }
1659  else
1660  {
1661  fnd()++;
1662  }
1663  }
1664  }
1665  }
1666 
1667  // Perform various checks on common points
1668 
1669  // Check all vertices shared (duplicate point)
1670  if (checkDuplicateFaces(faceI, nCommonPoints, nBaffleFaces, setPtr))
1671  {
1672  nErrorDuplicate++;
1673  }
1674 
1675  // Check common vertices are consecutive on both faces
1676  if (checkCommonOrder(faceI, nCommonPoints, setPtr))
1677  {
1678  nErrorOrder++;
1679  }
1680  }
1681 
1682  reduce(nBaffleFaces, sumOp<label>());
1683  reduce(nErrorDuplicate, sumOp<label>());
1684  reduce(nErrorOrder, sumOp<label>());
1685 
1686  if (nBaffleFaces)
1687  {
1688  Info<< " Number of identical duplicate faces (baffle faces): "
1689  << nBaffleFaces << endl;
1690  }
1691 
1692  if (nErrorDuplicate > 0 || nErrorOrder > 0)
1693  {
1694  // These are actually warnings, not errors.
1695  if (nErrorDuplicate > 0)
1696  {
1697  Info<< " <<Number of duplicate (not baffle) faces found: "
1698  << nErrorDuplicate
1699  << ". This might indicate a problem." << endl;
1700  }
1701 
1702  if (nErrorOrder > 0)
1703  {
1704  Info<< " <<Number of faces with non-consecutive shared points: "
1705  << nErrorOrder << ". This might indicate a problem." << endl;
1706  }
1707 
1708  return false; //return true;
1709  }
1710  else
1711  {
1712  if (debug || report)
1713  {
1714  Info<< " Face-face connectivity OK." << endl;
1715  }
1716 
1717  return false;
1718  }
1719 }
1720 
1721 
1722 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1723 
1724 bool Foam::primitiveMesh::checkClosedBoundary(const bool report) const
1725 {
1726  return checkClosedBoundary(faceAreas(), report, PackedBoolList(0));
1727 }
1728 
1729 
1732  const bool report,
1733  labelHashSet* setPtr,
1734  labelHashSet* aspectSetPtr,
1735  const Vector<label>& solutionD
1736 ) const
1737 {
1738  return checkClosedCells
1739  (
1740  faceAreas(),
1741  cellVolumes(),
1742  report,
1743  setPtr,
1744  aspectSetPtr,
1745  solutionD
1746  );
1747 }
1748 
1749 
1751 (
1752  const bool report,
1753  labelHashSet* setPtr
1754 ) const
1755 {
1756  return checkFaceAreas
1757  (
1758  faceAreas(),
1759  report,
1760  false, // detailedReport,
1761  setPtr
1762  );
1763 }
1764 
1765 
1767 (
1768  const bool report,
1769  labelHashSet* setPtr
1770 ) const
1771 {
1772  return checkCellVolumes
1773  (
1774  cellVolumes(),
1775  report,
1776  false, // detailedReport,
1777  setPtr
1778  );
1779 }
1780 
1781 
1783 (
1784  const bool report,
1785  labelHashSet* setPtr
1786 ) const
1787 {
1788  return checkFaceOrthogonality
1789  (
1790  faceAreas(),
1791  cellCentres(),
1792  report,
1793  setPtr
1794  );
1795 }
1796 
1797 
1799 (
1800  const bool report,
1801  const scalar minPyrVol,
1802  labelHashSet* setPtr
1803 ) const
1804 {
1805  return checkFacePyramids
1806  (
1807  points(),
1808  cellCentres(),
1809  report,
1810  false, // detailedReport,
1811  minPyrVol,
1812  setPtr
1813  );
1814 }
1815 
1816 
1818 (
1819  const bool report,
1820  labelHashSet* setPtr
1821 ) const
1822 {
1823  return checkFaceSkewness
1824  (
1825  points(),
1826  faceCentres(),
1827  faceAreas(),
1828  cellCentres(),
1829  report,
1830  setPtr
1831  );
1832 }
1833 
1834 
1837  const bool report,
1838  const scalar maxDeg,
1839  labelHashSet* setPtr
1840 ) const
1841 {
1842  return checkFaceAngles
1843  (
1844  points(),
1845  faceAreas(),
1846  report,
1847  maxDeg,
1848  setPtr
1849  );
1850 }
1851 
1852 
1855  const bool report,
1856  const scalar warnFlatness,
1857  labelHashSet* setPtr
1858 ) const
1859 {
1860  return checkFaceFlatness
1861  (
1862  points(),
1863  faceCentres(),
1864  faceAreas(),
1865  report,
1866  warnFlatness,
1867  setPtr
1868  );
1869 }
1870 
1871 
1874  const bool report,
1875  labelHashSet* setPtr
1876 ) const
1877 {
1878  return checkConcaveCells
1879  (
1880  faceAreas(),
1881  faceCentres(),
1882  report,
1883  setPtr
1884  );
1885 }
1886 
1887 
1888 bool Foam::primitiveMesh::checkTopology(const bool report) const
1889 {
1890  label noFailedChecks = 0;
1891 
1892  if (checkPoints(report)) noFailedChecks++;
1893  if (checkUpperTriangular(report)) noFailedChecks++;
1894  if (checkCellsZipUp(report)) noFailedChecks++;
1895  if (checkFaceVertices(report)) noFailedChecks++;
1896  if (checkFaceFaces(report)) noFailedChecks++;
1897 
1898  if (noFailedChecks == 0)
1899  {
1900  if (debug || report)
1901  {
1902  Info<< " Mesh topology OK." << endl;
1903  }
1904 
1905  return false;
1906  }
1907  else
1908  {
1909  if (debug || report)
1910  {
1911  Info<< " Failed " << noFailedChecks
1912  << " mesh topology checks." << endl;
1913  }
1914 
1915  return true;
1916  }
1917 }
1918 
1919 
1920 bool Foam::primitiveMesh::checkGeometry(const bool report) const
1921 {
1922  label noFailedChecks = 0;
1923 
1924  if (checkClosedBoundary(report)) noFailedChecks++;
1925  if (checkClosedCells(report)) noFailedChecks++;
1926  if (checkFaceAreas(report)) noFailedChecks++;
1927  if (checkCellVolumes(report)) noFailedChecks++;
1928  if (checkFaceOrthogonality(report)) noFailedChecks++;
1929  if (checkFacePyramids(report)) noFailedChecks++;
1930  if (checkFaceSkewness(report)) noFailedChecks++;
1931 
1932  if (noFailedChecks == 0)
1933  {
1934  if (debug || report)
1935  {
1936  Info<< " Mesh geometry OK." << endl;
1937  }
1938  return false;
1939  }
1940  else
1941  {
1942  if (debug || report)
1943  {
1944  Info<< " Failed " << noFailedChecks
1945  << " mesh geometry checks." << endl;
1946  }
1947 
1948  return true;
1949  }
1950 }
1951 
1952 
1953 bool Foam::primitiveMesh::checkMesh(const bool report) const
1954 {
1955  if (debug)
1956  {
1957  Info<< "bool primitiveMesh::checkMesh(const bool report) const: "
1958  << "checking primitiveMesh" << endl;
1959  }
1960 
1961  label noFailedChecks = checkTopology(report) + checkGeometry(report);
1962 
1963  if (noFailedChecks == 0)
1964  {
1965  if (debug || report)
1966  {
1967  Info<< "Mesh OK." << endl;
1968  }
1969 
1970  return false;
1971  }
1972  else
1973  {
1974  if (debug || report)
1975  {
1976  Info<< " Failed " << noFailedChecks
1977  << " mesh checks." << endl;
1978  }
1979 
1980  return true;
1981  }
1982 }
1983 
1984 
1985 Foam::scalar Foam::primitiveMesh::setClosedThreshold(const scalar val)
1986 {
1987  scalar prev = closedThreshold_;
1988  closedThreshold_ = val;
1989 
1990  return prev;
1991 }
1992 
1993 
1994 Foam::scalar Foam::primitiveMesh::setAspectThreshold(const scalar val)
1995 {
1996  scalar prev = aspectThreshold_;
1997  aspectThreshold_ = val;
1998 
1999  return prev;
2000 }
2001 
2002 
2003 Foam::scalar Foam::primitiveMesh::setNonOrthThreshold(const scalar val)
2004 {
2005  scalar prev = nonOrthThreshold_;
2006  nonOrthThreshold_ = val;
2007 
2008  return prev;
2009 }
2010 
2011 
2012 Foam::scalar Foam::primitiveMesh::setSkewThreshold(const scalar val)
2013 {
2014  scalar prev = skewThreshold_;
2015  skewThreshold_ = val;
2016 
2017  return prev;
2018 }
2019 
2020 
2021 // ************************************************************************* //
Foam::maxOp
Definition: ops.H:172
Foam::polyMeshGenChecks::checkClosedBoundary
bool checkClosedBoundary(const polyMeshGen &, const bool report=false)
Check boundary closedness.
Definition: polyMeshGenChecksGeometry.C:45
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::edgeList
List< edge > edgeList
Definition: edgeList.H:38
Foam::polyMeshGenChecks::checkFaceSkewness
void checkFaceSkewness(const polyMeshGen &, scalarField &, const boolList *changedFacePtr=NULL)
Check face skewness.
Definition: polyMeshGenChecksGeometry.C:1138
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::primitiveMesh::checkMesh
bool checkMesh(const bool report=false) const
Check mesh for correctness. Returns false for no error.
Definition: primitiveMeshCheck.C:2131
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
primitiveMeshTools.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
Foam::primitiveMesh::checkFaceVertices
bool checkFaceVertices(const bool report=false, labelHashSet *setPtr=NULL) const
Check uniqueness of face vertices.
Definition: primitiveMeshCheck.C:1553
Foam::labelMax
static const label labelMax
Definition: label.H:62
Foam::minOp
Definition: ops.H:173
Foam::primitiveMesh::faceOwner
virtual const labelList & faceOwner() const =0
Face face-owner addresing.
Foam::primitiveMesh::setSkewThreshold
static scalar setSkewThreshold(const scalar)
Set the skewness warning threshold as percentage.
Definition: primitiveMeshCheck.C:2012
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
Foam::primitiveMesh::checkDuplicateFaces
bool checkDuplicateFaces(const label, const Map< label > &, label &nBaffleFaces, labelHashSet *) const
Check if all points on face are shared with another face.
Definition: primitiveMeshCheck.C:1629
Foam::polyMeshGenChecks::checkUpperTriangular
bool checkUpperTriangular(const polyMeshGen &, const bool report=false, labelHashSet *setPtr=NULL)
Check face ordering.
Definition: polyMeshGenChecksTopology.C:118
Foam::checkGeometry
label checkGeometry(const polyMesh &mesh, const bool allGeometry, const autoPtr< surfaceWriter > &)
Definition: checkGeometry.C:478
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
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 > &)
Foam::HashSet< label, Hash< label > >
Foam::polyMeshGenChecks::checkFaceFlatness
bool checkFaceFlatness(const polyMeshGen &, const bool report, const scalar warnFlatness, labelHashSet *setPtr=NULL, const boolList *changedFacePtr=NULL)
Check face warpage: decompose face and check ratio between.
Definition: polyMeshGenChecksGeometry.C:1777
Foam::primitiveMesh::checkFaceAreas
bool checkFaceAreas(const bool report=false, labelHashSet *setPtr=NULL) const
Check for negative face areas.
Definition: primitiveMeshCheck.C:338
Foam::polyMeshGenChecks::checkFaceVertices
bool checkFaceVertices(const polyMeshGen &, const bool report=false, labelHashSet *setPtr=NULL)
Check uniqueness of face vertices.
Definition: polyMeshGenChecksTopology.C:373
Foam::primitiveMesh::checkTopology
bool checkTopology(const bool report=false) const
Check mesh topology for correctness.
Definition: primitiveMeshCheck.C:2065
Foam::primitiveMesh::faceNeighbour
virtual const labelList & faceNeighbour() const =0
Face face-neighbour addressing.
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::primitiveMesh::checkConcaveCells
bool checkConcaveCells(const vectorField &fAreas, const pointField &fCentres, const bool report, labelHashSet *setPtr) const
Check for concave cells by the planes of faces.
Definition: primitiveMeshCheck.C:874
Foam::polyMeshGenChecks::checkFaceAngles
bool checkFaceAngles(const polyMeshGen &, const bool report=false, const scalar maxDeg=10, labelHashSet *setPtr=NULL, const boolList *changedFacePtr=NULL)
Check face angles.
Definition: polyMeshGenChecksGeometry.C:1620
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:359
Foam::primitiveMesh::setClosedThreshold
static scalar setClosedThreshold(const scalar)
Set the closedness ratio warning threshold.
Definition: primitiveMeshCheck.C:1985
Foam::cmptMax
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:224
SortableList.H
Foam::primitiveMesh::checkGeometry
bool checkGeometry(const bool report=false) const
Check mesh geometry (& implicitly topology) for correctness.
Definition: primitiveMeshCheck.C:2098
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::primitiveMesh::setNonOrthThreshold
static scalar setNonOrthThreshold(const scalar)
Set the non-orthogonality warning threshold in degrees.
Definition: primitiveMeshCheck.C:2003
Foam::primitiveMesh::checkFaceFlatness
bool checkFaceFlatness(const bool report, labelHashSet *setPtr) const
Check face warpage: decompose face and check ratio between.
Definition: primitiveMeshCheck.C:1020
Foam::Map< label >::iterator
HashTable< T, label, Hash< label > >::iterator iterator
Definition: Map.H:56
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::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyMeshGenChecks::checkCellsZipUp
bool checkCellsZipUp(const polyMeshGen &, const bool report=false, labelHashSet *setPtr=NULL)
Check cell zip-up.
Definition: polyMeshGenChecksTopology.C:254
Foam::cellList
List< cell > cellList
list of cells
Definition: cellList.H:42
Foam::primitiveMesh::checkFaceOrthogonality
bool checkFaceOrthogonality(const bool report=false, labelHashSet *setPtr=NULL) const
Check for non-orthogonality.
Definition: primitiveMeshCheck.C:464
Foam::polyMeshGenChecks::checkPoints
bool checkPoints(const polyMeshGen &, const bool report=false, labelHashSet *setPtr=NULL)
Check for unused points.
Definition: polyMeshGenChecksTopology.C:46
Foam::FatalError
error FatalError
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
Foam::primitiveMesh::aspectThreshold_
static const debug::tolerancesSwitch aspectThreshold_
Aspect ratio warning threshold.
Definition: primitiveMesh.H:327
Foam::primitiveMesh::checkClosedBoundary
bool checkClosedBoundary(const bool report=false) const
Check boundary for closedness.
Definition: primitiveMeshCheck.C:79
Foam::primitiveMesh::checkFaceSkewness
bool checkFaceSkewness(const bool report=false, labelHashSet *setPtr=NULL) const
Check face skewness.
Definition: primitiveMeshCheck.C:669
Foam::primitiveMesh::checkUpperTriangular
bool checkUpperTriangular(const bool report=false, labelHashSet *setPtr=NULL) const
Check face ordering.
Definition: primitiveMeshCheck.C:1285
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::primitiveMesh::closedThreshold_
static const debug::tolerancesSwitch closedThreshold_
Static data to control mesh checking.
Definition: primitiveMesh.H:324
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
Foam::primitiveMesh::checkClosedCells
bool checkClosedCells(const bool report=false, labelHashSet *setPtr=NULL, labelHashSet *highAspectSetPtr=NULL) const
Check cells for closedness.
Definition: primitiveMeshCheck.C:135
Foam::primitiveMesh::nonOrthThreshold_
static Foam::debug::tolerancesSwitch nonOrthThreshold_
Non-orthogonality warning threshold in deg.
Definition: primitiveMesh.H:330
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
Foam::primitiveMesh::checkFaceAngles
bool checkFaceAngles(const bool report=false, labelHashSet *setPtr=NULL) const
Check face angles.
Definition: primitiveMeshCheck.C:892
Foam::primitiveMesh::planarCosAngle_
static scalar planarCosAngle_
Threshold where faces are considered coplanar.
Definition: primitiveMesh.H:246
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::primitiveMesh::setAspectThreshold
static scalar setAspectThreshold(const scalar)
Set the aspect ratio warning threshold.
Definition: primitiveMeshCheck.C:1994
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::primitiveMesh::skewThreshold_
static const debug::tolerancesSwitch skewThreshold_
Skewness warning threshold.
Definition: primitiveMesh.H:333
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
Foam::Vector< scalar >
Foam::primitiveMesh::checkCellVolumes
bool checkCellVolumes(const bool report=false, labelHashSet *setPtr=NULL) const
Check for negative cell volumes.
Definition: primitiveMeshCheck.C:397
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::primitiveMesh::checkPoints
bool checkPoints(const bool report=false, labelHashSet *setPtr=NULL) const
Check for unused points.
Definition: primitiveMeshCheck.C:814
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:259
Foam::primitiveMesh::checkFaceFaces
bool checkFaceFaces(const bool report=false, labelHashSet *setPtr=NULL) const
Check face-face connectivity.
Definition: primitiveMeshCheck.C:1834
points
const pointField & points
Definition: gmvOutputHeader.H:1
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
EdgeMap.H
Foam::radToDeg
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
Definition: unitConversion.H:51
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
ListOps.H
Various functions to operate on Lists.
Foam::PackedList::size
label size() const
Number of entries.
Definition: PackedListI.H:714
Foam::primitiveMesh::checkCellsZipUp
bool checkCellsZipUp(const bool report=false, labelHashSet *setPtr=NULL) const
Check cell zip-up.
Definition: primitiveMeshCheck.C:1452
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
Foam::primitiveMesh::checkFacePyramids
bool checkFacePyramids(const bool report=false, const scalar minPyrVol=-SMALL, labelHashSet *setPtr=NULL) const
Check face pyramid volume.
Definition: primitiveMeshCheck.C:585
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
pyramidPointFaceRef.H
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::polyMeshGenChecks::checkClosedCells
bool checkClosedCells(const polyMeshGen &, const bool report=false, const scalar aspectWarn=1000, labelHashSet *setPtr=NULL)
Check cells for closedness.
Definition: polyMeshGenChecksGeometry.C:113
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:258
Foam::checkTopology
label checkTopology(const polyMesh &, const bool, const bool, const autoPtr< surfaceWriter > &)
Definition: checkTopology.C:42
Foam::primitiveMesh::checkCommonOrder
bool checkCommonOrder(const label, const Map< label > &, labelHashSet *) const
Check that shared points are in consecutive order.
Definition: primitiveMeshCheck.C:1671
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