checkGeometry.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 "meshCheck.H"
27 #include "PatchTools.H"
28 #include "polyMeshCheck.H"
29 #include "cellSet.H"
30 #include "faceSet.H"
31 #include "pointSet.H"
32 #include "EdgeMap.H"
33 #include "wedgePolyPatch.H"
35 
36 #include "vtkSurfaceWriter.H"
37 #include "setWriter.H"
38 #include "writeFile.H"
40 
41 #include "mergeAndWrite.H"
42 #include "Time.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
47 (
48  const polyMesh& mesh,
49  const wedgePolyPatch& wpp
50 )
51 {
53 
54  scalar wppCosAngle = wpp.cosAngle();
55 
57  {
58  if
59  (
60  patchi != wpp.index()
61  && patches[patchi].size()
62  && isA<wedgePolyPatch>(patches[patchi])
63  )
64  {
65  const wedgePolyPatch& pp =
66  refCast<const wedgePolyPatch>(patches[patchi]);
67 
68  // Calculate (cos of) angle to wpp (not pp!) centre normal
69  scalar ppCosAngle = wpp.centreNormal() & pp.n();
70 
71  if
72  (
73  pp.size() == wpp.size()
74  && mag(pp.axis() & wpp.axis()) >= (1-1e-3)
75  && mag(ppCosAngle - wppCosAngle) >= 1e-3
76  )
77  {
78  return patchi;
79  }
80  }
81  }
82  return -1;
83 }
84 
85 
87 (
88  const polyMesh& mesh,
89  const bool report,
91  labelHashSet* setPtr
92 )
93 {
94  // To mark edges without calculating edge addressing
95  EdgeMap<label> edgesInError;
96 
97  const pointField& p = mesh.points();
98  const faceList& fcs = mesh.faces();
99 
100 
103  {
104  if (patches[patchi].size() && isA<wedgePolyPatch>(patches[patchi]))
105  {
106  const wedgePolyPatch& pp =
107  refCast<const wedgePolyPatch>(patches[patchi]);
108 
109  scalar wedgeAngle = acos(pp.cosAngle());
110 
111  if (report)
112  {
113  Info<< " Wedge " << pp.name() << " with angle "
114  << radToDeg(wedgeAngle) << " degrees"
115  << endl;
116  }
117 
118  // Find opposite
119  label oppositePatchi = findOppositeWedge(mesh, pp);
120 
121  if (oppositePatchi == -1)
122  {
123  if (report)
124  {
125  Info<< " ***Cannot find opposite wedge for wedge "
126  << pp.name() << endl;
127  }
128  return true;
129  }
130 
131  const wedgePolyPatch& opp =
132  refCast<const wedgePolyPatch>(patches[oppositePatchi]);
133 
134 
135  if (mag(opp.axis() & pp.axis()) < (1-1e-3))
136  {
137  if (report)
138  {
139  Info<< " ***Wedges do not have the same axis."
140  << " Encountered " << pp.axis()
141  << " on patch " << pp.name()
142  << " which differs from " << opp.axis()
143  << " on opposite wedge patch" << opp.axis()
144  << endl;
145  }
146  return true;
147  }
148 
149 
150 
151  // Mark edges on wedgePatches
152  forAll(pp, i)
153  {
154  const face& f = pp[i];
155  forAll(f, fp)
156  {
157  label p0 = f[fp];
158  label p1 = f.nextLabel(fp);
159  edgesInError.insert(edge(p0, p1), -1); // non-error value
160  }
161  }
162 
163 
164  // Check that wedge patch is flat
165  const point& p0 = p[pp.meshPoints()[0]];
166  forAll(pp.meshPoints(), i)
167  {
168  const point& pt = p[pp.meshPoints()[i]];
169  scalar d = mag((pt - p0) & pp.n());
170 
171  if (d > rootSmall)
172  {
173  if (report)
174  {
175  Info<< " ***Wedge patch " << pp.name() << " not planar."
176  << " Point " << pt << " is not in patch plane by "
177  << d << " metre."
178  << endl;
179  }
180  return true;
181  }
182  }
183  }
184  }
185 
186 
187 
188  // Check all non-wedge faces
189  label nEdgesInError = 0;
190 
191  forAll(fcs, facei)
192  {
193  const face& f = fcs[facei];
194 
195  forAll(f, fp)
196  {
197  label p0 = f[fp];
198  label p1 = f.nextLabel(fp);
199  if (p0 < p1)
200  {
201  vector d(p[p1]-p[p0]);
202  scalar magD = mag(d);
203 
204  if (magD > rootVSmall)
205  {
206  d /= magD;
207 
208  // Check how many empty directions are used by the edge.
209  label nEmptyDirs = 0;
210  label nNonEmptyDirs = 0;
211  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
212  {
213  if (mag(d[cmpt]) > 1e-6)
214  {
215  if (directions[cmpt] == 0)
216  {
217  nEmptyDirs++;
218  }
219  else
220  {
221  nNonEmptyDirs++;
222  }
223  }
224  }
225 
226  if (nEmptyDirs == 0)
227  {
228  // Purely in ok directions.
229  }
230  else if (nEmptyDirs == 1)
231  {
232  // Ok if purely in empty directions.
233  if (nNonEmptyDirs > 0)
234  {
235  if (edgesInError.insert(edge(p0, p1), facei))
236  {
237  nEdgesInError++;
238  }
239  }
240  }
241  else if (nEmptyDirs > 1)
242  {
243  // Always an error
244  if (edgesInError.insert(edge(p0, p1), facei))
245  {
246  nEdgesInError++;
247  }
248  }
249  }
250  }
251  }
252  }
253 
254  label nErrorEdges = returnReduce(nEdgesInError, sumOp<label>());
255 
256  if (nErrorEdges > 0)
257  {
258  if (report)
259  {
260  Info<< " ***Number of edges not aligned with or perpendicular to "
261  << "non-empty directions: " << nErrorEdges << endl;
262  }
263 
264  if (setPtr)
265  {
266  setPtr->resize(2*nEdgesInError);
267  forAllConstIter(EdgeMap<label>, edgesInError, iter)
268  {
269  if (iter() >= 0)
270  {
271  setPtr->insert(iter.key()[0]);
272  setPtr->insert(iter.key()[1]);
273  }
274  }
275  }
276 
277  return true;
278  }
279  else
280  {
281  if (report)
282  {
283  Info<< " All edges aligned with or perpendicular to "
284  << "non-empty directions." << endl;
285  }
286  return false;
287  }
288 }
289 
290 
291 namespace Foam
292 {
293  //- Default transformation behaviour for position
295  {
296  public:
297 
298  //- Transform patch-based field
299  void operator()
300  (
301  const coupledPolyPatch& cpp,
302  List<pointField>& pts
303  ) const
304  {
305  // Each element of pts is all the points in the face. Convert into
306  // lists of size cpp to transform.
307 
308  List<pointField> newPts(pts.size());
309  forAll(pts, facei)
310  {
311  newPts[facei].setSize(pts[facei].size());
312  }
313 
314  label index = 0;
315  while (true)
316  {
317  label n = 0;
318 
319  // Extract for every face the i'th position
320  pointField ptsAtIndex(pts.size(), Zero);
321  forAll(cpp, facei)
322  {
323  const pointField& facePts = pts[facei];
324  if (facePts.size() > index)
325  {
326  ptsAtIndex[facei] = facePts[index];
327  n++;
328  }
329  }
330 
331  if (n == 0)
332  {
333  break;
334  }
335 
336  // Now ptsAtIndex will have for every face either zero or
337  // the position of the i'th vertex. Transform.
338  cpp.transform().transformPosition(ptsAtIndex, ptsAtIndex);
339 
340  // Extract back from ptsAtIndex into newPts
341  forAll(cpp, facei)
342  {
343  pointField& facePts = newPts[facei];
344  if (facePts.size() > index)
345  {
346  facePts[index] = ptsAtIndex[facei];
347  }
348  }
349 
350  index++;
351  }
352 
353  pts.transfer(newPts);
354  }
355  };
356 }
357 
358 
360 (
361  const polyMesh& mesh,
362  const bool report,
363  labelHashSet* setPtr
364 )
365 {
366  const pointField& p = mesh.points();
367  const faceList& fcs = mesh.faces();
369 
370  // Zero'th point on coupled faces
371  // pointField nbrZeroPoint(fcs.size()-mesh.nInternalFaces(), vector::max);
372  List<pointField> nbrPoints(fcs.size() - mesh.nInternalFaces());
373 
374  // Exchange zero point
376  {
377  if (patches[patchi].coupled())
378  {
379  const coupledPolyPatch& cpp = refCast<const coupledPolyPatch>
380  (
381  patches[patchi]
382  );
383 
384  forAll(cpp, i)
385  {
386  label bFacei = cpp.start() + i - mesh.nInternalFaces();
387  const face& f = cpp[i];
388  nbrPoints[bFacei].setSize(f.size());
389  forAll(f, fp)
390  {
391  const point& p0 = p[f[fp]];
392  nbrPoints[bFacei][fp] = p0;
393  }
394  }
395  }
396  }
397  syncTools::syncBoundaryFaceList
398  (
399  mesh,
400  nbrPoints,
403  );
404 
405  // Compare to local ones. Use same tolerance as for matching
406  label nErrorFaces = 0;
407  scalar avgMismatch = 0;
408  label nCoupledPoints = 0;
409 
411  {
412  if (patches[patchi].coupled())
413  {
414  const coupledPolyPatch& cpp =
415  refCast<const coupledPolyPatch>(patches[patchi]);
416 
417  if (cpp.owner())
418  {
419  scalarField smallDist
420  (
421  cpp.calcFaceTol
422  (
423  // cpp.matchTolerance(),
424  cpp,
425  cpp.points(),
426  cpp.faceCentres()
427  )
428  );
429 
430  forAll(cpp, i)
431  {
432  label bFacei = cpp.start() + i - mesh.nInternalFaces();
433  const face& f = cpp[i];
434 
435  if (f.size() != nbrPoints[bFacei].size())
436  {
438  << "Local face size : " << f.size()
439  << " does not equal neighbour face size : "
440  << nbrPoints[bFacei].size()
441  << abort(FatalError);
442  }
443 
444  label fp = 0;
445  forAll(f, j)
446  {
447  const point& p0 = p[f[fp]];
448  scalar d = mag(p0 - nbrPoints[bFacei][j]);
449 
450  if (d > smallDist[i])
451  {
452  if (setPtr)
453  {
454  setPtr->insert(cpp.start()+i);
455  }
456  nErrorFaces++;
457 
458  break;
459  }
460 
461  avgMismatch += d;
462  nCoupledPoints++;
463 
464  fp = f.rcIndex(fp);
465  }
466  }
467  }
468  }
469  }
470 
471  reduce(nErrorFaces, sumOp<label>());
472  reduce(avgMismatch, maxOp<scalar>());
473  reduce(nCoupledPoints, sumOp<label>());
474 
475  if (nCoupledPoints > 0)
476  {
477  avgMismatch /= nCoupledPoints;
478  }
479 
480  if (nErrorFaces > 0)
481  {
482  if (report)
483  {
484  Info<< " **Error in coupled point location: "
485  << nErrorFaces
486  << " faces have their 0th or consecutive vertex not opposite"
487  << " their coupled equivalent. Average mismatch "
488  << avgMismatch << "."
489  << endl;
490  }
491 
492  return true;
493  }
494  else
495  {
496  if (report)
497  {
498  Info<< " Coupled point location match (average "
499  << avgMismatch << ") OK." << endl;
500  }
501 
502  return false;
503  }
504 }
505 
506 
508 (
509  const polyMesh& mesh,
510  const bool allGeometry,
511  const scalar nonOrthThreshold,
512  const scalar skewThreshold,
513  const autoPtr<surfaceWriter>& surfWriter,
515 )
516 {
517  const scalar closedThreshold = 1.0e-6;
518  const scalar aspectThreshold = 1000;
519  const scalar planarCosAngle = 1.0e-6;
520 
521  label noFailedChecks = 0;
522 
523  Info<< "\nChecking geometry..." << endl;
524 
525  // Get a small relative length from the bounding box
526  const boundBox& globalBb = mesh.bounds();
527 
528  Info<< " Overall domain bounding box "
529  << globalBb.min() << " " << globalBb.max() << endl;
530 
531  if (allGeometry)
532  {
533  Info<< " Patch bounding boxes" << endl;
534 
536 
538  {
539  const polyPatch& pp = patches[patchi];
540 
541  if (!isA<processorPolyPatch>(pp))
542  {
543  Info<< " " << setw(20) << pp.name();
544 
545  const pointField& pts = pp.points();
546  const labelList& mp = pp.meshPoints();
547 
548  if (returnReduce(mp.size(), sumOp<label>()) > 0)
549  {
551  forAll(mp, i)
552  {
553  bb.min() = min(bb.min(), pts[mp[i]]);
554  bb.max() = max(bb.max(), pts[mp[i]]);
555  }
556  reduce(bb.min(), minOp<vector>());
557  reduce(bb.max(), maxOp<vector>());
558  Info<< ' ' << bb;
559  }
560  }
561  Info<< endl;
562  }
563  }
564 
565  // Min length
566  scalar minDistSqr = magSqr(1e-6 * globalBb.span());
567 
568  // Geometric directions
569  const Vector<label> validDirs = (mesh.geometricD() + Vector<label>::one)/2;
570  Info<< " Mesh has " << mesh.nGeometricD()
571  << " geometric (non-empty/wedge) directions " << validDirs << endl;
572 
573  // Solution directions
574  const Vector<label> solDirs = (mesh.solutionD() + Vector<label>::one)/2;
575  Info<< " Mesh has " << mesh.nSolutionD()
576  << " solution (non-empty) directions " << solDirs << endl;
577 
578  if (mesh.nGeometricD() < 3)
579  {
580  pointSet nonAlignedPoints(mesh, "nonAlignedEdges", mesh.nPoints()/100);
581 
582  if
583  (
584  (
585  validDirs != solDirs
586  && checkWedges(mesh, true, validDirs, &nonAlignedPoints)
587  )
588  || (
589  validDirs == solDirs
591  (
592  mesh,
593  true,
594  validDirs,
595  &nonAlignedPoints
596  )
597  )
598  )
599  {
600  noFailedChecks++;
601  label nNonAligned = returnReduce
602  (
603  nonAlignedPoints.size(),
604  sumOp<label>()
605  );
606 
607  if (nNonAligned > 0)
608  {
609  if (setWriter.valid())
610  {
611  Info<< " <<Writing " << nNonAligned
612  << " points on non-aligned edges to set "
613  << nonAlignedPoints.name() << endl;
614  nonAlignedPoints.instance() = mesh.pointsInstance();
615  nonAlignedPoints.write();
616  meshCheck::mergeAndWrite(setWriter, nonAlignedPoints);
617  }
618  }
619  }
620  }
621 
622  if (meshCheck::checkClosedBoundary(mesh, true)) noFailedChecks++;
623 
624  {
625  cellSet cells(mesh, "nonClosedCells", mesh.nCells()/100+1);
626  cellSet aspectCells(mesh, "highAspectRatioCells", mesh.nCells()/100+1);
627  if
628  (
630  (
631  mesh,
632  closedThreshold,
633  aspectThreshold,
634  true,
635  &cells,
636  &aspectCells,
637  mesh.geometricD()
638  )
639  )
640  {
641  noFailedChecks++;
642 
643  label nNonClosed = returnReduce(cells.size(), sumOp<label>());
644 
645  if (nNonClosed > 0)
646  {
647  if (setWriter.valid())
648  {
649  Info<< " <<Writing " << nNonClosed
650  << " non closed cells to set " << cells.name() << endl;
651  cells.instance() = mesh.pointsInstance();
652  cells.write();
653  }
654 
655  if (surfWriter.valid())
656  {
657  meshCheck::mergeAndWrite(surfWriter(), cells);
658  }
659  }
660  }
661 
662  label nHighAspect = returnReduce(aspectCells.size(), sumOp<label>());
663 
664  if (nHighAspect > 0)
665  {
666  if (setWriter.valid())
667  {
668  Info<< " <<Writing " << nHighAspect
669  << " cells with high aspect ratio to set "
670  << aspectCells.name() << endl;
671  aspectCells.instance() = mesh.pointsInstance();
672  aspectCells.write();
673  }
674 
675  if (surfWriter.valid())
676  {
677  meshCheck::mergeAndWrite(surfWriter(), aspectCells);
678  }
679  }
680  }
681 
682  {
683  faceSet faces(mesh, "zeroAreaFaces", mesh.nFaces()/100+1);
684  if (meshCheck::checkFaceAreas(mesh, true, &faces))
685  {
686  noFailedChecks++;
687 
688  label nFaces = returnReduce(faces.size(), sumOp<label>());
689 
690  if (nFaces > 0)
691  {
692  if (setWriter.valid())
693  {
694  Info<< " <<Writing " << nFaces
695  << " zero area faces to set " << faces.name() << endl;
696  faces.instance() = mesh.pointsInstance();
697  faces.write();
698  }
699 
700  if (surfWriter.valid())
701  {
702  meshCheck::mergeAndWrite(surfWriter(), faces);
703  }
704  }
705  }
706  }
707 
708  {
709  cellSet cells(mesh, "zeroVolumeCells", mesh.nCells()/100+1);
711  {
712  noFailedChecks++;
713 
714  label nCells = returnReduce(cells.size(), sumOp<label>());
715 
716  if (nCells > 0)
717  {
718  if (setWriter.valid())
719  {
720  Info<< " <<Writing " << nCells
721  << " zero volume cells to set " << cells.name() << endl;
722  cells.instance() = mesh.pointsInstance();
723  cells.write();
724  }
725 
726  if (surfWriter.valid())
727  {
728  meshCheck::mergeAndWrite(surfWriter(), cells);
729  }
730  }
731  }
732  }
733 
734  {
735  faceSet faces(mesh, "nonOrthoFaces", mesh.nFaces()/100+1);
736  if
737  (
739  (
740  mesh,
741  nonOrthThreshold,
742  true,
743  &faces
744  )
745  )
746  {
747  noFailedChecks++;
748  }
749 
750  label nFaces = returnReduce(faces.size(), sumOp<label>());
751 
752  if (nFaces > 0)
753  {
754  if (setWriter.valid())
755  {
756  Info<< " <<Writing " << nFaces
757  << " non-orthogonal faces to set " << faces.name() << endl;
758  faces.instance() = mesh.pointsInstance();
759  faces.write();
760  }
761 
762  if (surfWriter.valid())
763  {
764  meshCheck::mergeAndWrite(surfWriter(), faces);
765  }
766  }
767  }
768 
769  {
770  faceSet faces(mesh, "wrongOrientedFaces", mesh.nFaces()/100 + 1);
771  if (meshCheck::checkFacePyramids(mesh, true, -small, &faces))
772  {
773  noFailedChecks++;
774 
775  label nFaces = returnReduce(faces.size(), sumOp<label>());
776 
777  if (nFaces > 0)
778  {
779  if (setWriter.valid())
780  {
781  Info<< " <<Writing " << nFaces
782  << " faces with incorrect orientation to set "
783  << faces.name() << endl;
784  faces.instance() = mesh.pointsInstance();
785  faces.write();
786  }
787 
788  if (surfWriter.valid())
789  {
790  meshCheck::mergeAndWrite(surfWriter(), faces);
791  }
792  }
793  }
794  }
795 
796  {
797  faceSet faces(mesh, "skewFaces", mesh.nFaces()/100+1);
798  if
799  (
801  (
802  mesh,
803  skewThreshold,
804  true,
805  &faces
806  )
807  )
808  {
809  noFailedChecks++;
810 
811  label nFaces = returnReduce(faces.size(), sumOp<label>());
812 
813  if (nFaces > 0)
814  {
815  if (setWriter.valid())
816  {
817  Info<< " <<Writing " << nFaces
818  << " skew faces to set " << faces.name() << endl;
819  faces.instance() = mesh.pointsInstance();
820  faces.write();
821  }
822 
823  if (surfWriter.valid())
824  {
825  meshCheck::mergeAndWrite(surfWriter(), faces);
826  }
827  }
828  }
829  }
830 
831  {
832  faceSet faces(mesh, "coupledFaces", mesh.nFaces()/100 + 1);
833  if (checkCoupledPoints(mesh, true, &faces))
834  {
835  noFailedChecks++;
836 
837  label nFaces = returnReduce(faces.size(), sumOp<label>());
838 
839  if (nFaces > 0)
840  {
841  if (setWriter.valid())
842  {
843  Info<< " <<Writing " << nFaces
844  << " faces with incorrectly matched 0th "
845  "(or consecutive) vertex to set "
846  << faces.name() << endl;
847  faces.instance() = mesh.pointsInstance();
848  faces.write();
849  }
850 
851  if (surfWriter.valid())
852  {
853  meshCheck::mergeAndWrite(surfWriter(), faces);
854  }
855  }
856  }
857  }
858 
859  if (allGeometry)
860  {
861  faceSet faces(mesh, "lowQualityTetFaces", mesh.nFaces()/100+1);
862  if
863  (
865  (
866  mesh,
867  polyMeshTetDecomposition::minTetQuality,
868  true,
869  &faces
870  )
871  )
872  {
873  noFailedChecks++;
874 
875  label nFaces = returnReduce(faces.size(), sumOp<label>());
876 
877  if (nFaces > 0)
878  {
879  if (setWriter.valid())
880  {
881  Info<< " <<Writing " << nFaces
882  << " faces with low quality or negative volume "
883  << "decomposition tets to set " << faces.name() << endl;
884  faces.instance() = mesh.pointsInstance();
885  faces.write();
886  }
887 
888  if (surfWriter.valid())
889  {
890  meshCheck::mergeAndWrite(surfWriter(), faces);
891  }
892  }
893  }
894  }
895 
896  if (allGeometry)
897  {
898  // Note use of nPoints since don't want edge construction.
899  pointSet points(mesh, "shortEdges", mesh.nPoints()/1000 + 1);
900  if (meshCheck::checkEdgeLength(mesh, true, minDistSqr, &points))
901  {
902  // noFailedChecks++;
903 
905 
906  if (nPoints > 0)
907  {
908  if (setWriter.valid())
909  {
910  Info<< " <<Writing " << nPoints
911  << " points on short edges to set " << points.name()
912  << endl;
913  points.instance() = mesh.pointsInstance();
914  points.write();
916  }
917  }
918  }
919 
920  label nEdgeClose = returnReduce(points.size(), sumOp<label>());
921 
922  if (meshCheck::checkPointNearness(mesh, false, minDistSqr, &points))
923  {
924  // noFailedChecks++;
925 
927 
928  if (nPoints > nEdgeClose)
929  {
930  if (setWriter.valid())
931  {
932  pointSet nearPoints(mesh, "nearPoints", points);
933  Info<< " <<Writing " << nPoints
934  << " near (closer than " << Foam::sqrt(minDistSqr)
935  << " apart) points to set " << nearPoints.name()
936  << endl;
937  nearPoints.instance() = mesh.pointsInstance();
938  nearPoints.write();
939  meshCheck::mergeAndWrite(setWriter, nearPoints);
940  }
941  }
942  }
943  }
944 
945  if (allGeometry)
946  {
947  faceSet faces(mesh, "concaveFaces", mesh.nFaces()/100 + 1);
948  if (meshCheck::checkFaceAngles(mesh, true, degToRad(10), &faces))
949  {
950  // noFailedChecks++;
951 
952  label nFaces = returnReduce(faces.size(), sumOp<label>());
953 
954  if (nFaces > 0)
955  {
956  if (setWriter.valid())
957  {
958  Info<< " <<Writing " << nFaces
959  << " faces with concave angles to set " << faces.name()
960  << endl;
961  faces.instance() = mesh.pointsInstance();
962  faces.write();
963  }
964 
965  if (surfWriter.valid())
966  {
967  meshCheck::mergeAndWrite(surfWriter(), faces);
968  }
969  }
970  }
971  }
972 
973  if (allGeometry)
974  {
975  faceSet faces(mesh, "warpedFaces", mesh.nFaces()/100 + 1);
976  if (meshCheck::checkFaceFlatness(mesh, true, 0.8, &faces))
977  {
978  // noFailedChecks++;
979 
980  label nFaces = returnReduce(faces.size(), sumOp<label>());
981 
982  if (nFaces > 0)
983  {
984  if (setWriter.valid())
985  {
986  Info<< " <<Writing " << nFaces
987  << " warped faces to set " << faces.name() << endl;
988  faces.instance() = mesh.pointsInstance();
989  faces.write();
990  }
991 
992  if (surfWriter.valid())
993  {
994  meshCheck::mergeAndWrite(surfWriter(), faces);
995  }
996  }
997  }
998  }
999 
1000  if (allGeometry)
1001  {
1002  cellSet cells(mesh, "underdeterminedCells", mesh.nCells()/100);
1004  {
1005  noFailedChecks++;
1006 
1007  if (setWriter.valid())
1008  {
1009  label nCells = returnReduce(cells.size(), sumOp<label>());
1010 
1011  Info<< " <<Writing " << nCells
1012  << " under-determined cells to set " << cells.name()
1013  << endl;
1014  cells.instance() = mesh.pointsInstance();
1015  cells.write();
1016  }
1017 
1018  if (surfWriter.valid())
1019  {
1020  meshCheck::mergeAndWrite(surfWriter(), cells);
1021  }
1022  }
1023  }
1024 
1025  if (allGeometry)
1026  {
1027  cellSet cells(mesh, "concaveCells", mesh.nCells()/100);
1028  if (meshCheck::checkConcaveCells(mesh, planarCosAngle, true, &cells))
1029  {
1030  noFailedChecks++;
1031 
1032  if (setWriter.valid())
1033  {
1034  label nCells = returnReduce(cells.size(), sumOp<label>());
1035 
1036  Info<< " <<Writing " << nCells
1037  << " concave cells to set " << cells.name() << endl;
1038  cells.instance() = mesh.pointsInstance();
1039  cells.write();
1040  }
1041 
1042  if (surfWriter.valid())
1043  {
1044  meshCheck::mergeAndWrite(surfWriter(), cells);
1045  }
1046  }
1047  }
1048 
1049  if (allGeometry)
1050  {
1051  faceSet faces(mesh, "lowWeightFaces", mesh.nFaces()/100);
1052  if (meshCheck::checkFaceWeight(mesh, true, 0.05, &faces))
1053  {
1054  noFailedChecks++;
1055 
1056  if (setWriter.valid())
1057  {
1058  label nFaces = returnReduce(faces.size(), sumOp<label>());
1059 
1060  Info<< " <<Writing " << nFaces
1061  << " faces with low interpolation weights to set "
1062  << faces.name() << endl;
1063  faces.instance() = mesh.pointsInstance();
1064  faces.write();
1065  }
1066 
1067  if (surfWriter.valid())
1068  {
1069  meshCheck::mergeAndWrite(surfWriter(), faces);
1070  }
1071  }
1072  }
1073 
1074  if (allGeometry)
1075  {
1076  faceSet faces(mesh, "lowVolRatioFaces", mesh.nFaces()/100);
1077  if (meshCheck::checkVolRatio(mesh, true, 0.01, &faces))
1078  {
1079  noFailedChecks++;
1080 
1081  if (setWriter.valid())
1082  {
1083  label nFaces = returnReduce(faces.size(), sumOp<label>());
1084 
1085  Info<< " <<Writing " << nFaces
1086  << " faces with low volume ratio cells to set "
1087  << faces.name() << endl;
1088  faces.instance() = mesh.pointsInstance();
1089  faces.write();
1090  }
1091 
1092  if (surfWriter.valid())
1093  {
1094  meshCheck::mergeAndWrite(surfWriter(), faces);
1095  }
1096  }
1097  }
1098 
1099  if (allGeometry)
1100  {
1101  const fileName outputPath =
1102  mesh.time().globalPath()
1103  /functionObjects::writeFile::outputPrefix
1104  /(mesh.name() != polyMesh::defaultRegion ? mesh.name() : word())
1105  /"checkMesh"
1106  /mesh.time().name();
1107 
1109 
1110  // Compute coverage for all orig patches
1111  PtrList<scalarField> patchCoverage(patches.size());
1112  forAll(patches, nccPatchi)
1113  {
1114  if (isA<nonConformalCyclicPolyPatch>(patches[nccPatchi]))
1115  {
1116  const nonConformalCyclicPolyPatch& nccPp =
1117  refCast<const nonConformalCyclicPolyPatch>
1118  (patches[nccPatchi]);
1119 
1120  if (nccPp.owner())
1121  {
1122  const polyPatch& origPp = nccPp.origPatch();
1123  const polyPatch& nbrOrigPp = nccPp.nbrPatch().origPatch();
1124 
1126  nccPp.intersection();
1127 
1128  if (!patchCoverage.set(origPp.index()))
1129  {
1130  patchCoverage.set
1131  (
1132  origPp.index(),
1133  scalarField(origPp.size(), 0)
1134  );
1135  }
1136 
1137  patchCoverage[origPp.index()] +=
1138  intersection.srcCoverage();
1139 
1140  if (!patchCoverage.set(nbrOrigPp.index()))
1141  {
1142  patchCoverage.set
1143  (
1144  nbrOrigPp.index(),
1145  scalarField(nbrOrigPp.size(), 0)
1146  );
1147  }
1148 
1149  patchCoverage[nbrOrigPp.index()] +=
1150  intersection.tgtCoverage();
1151  }
1152  }
1153  }
1154 
1155  // Write out to surface files
1157  {
1158  if (patchCoverage.set(patchi))
1159  {
1160  const polyPatch& patch = patches[patchi];
1161 
1162  // Collect the patch geometry
1163  labelList pointToGlobal;
1164  labelList uniqueMeshPointLabels;
1166  autoPtr<globalIndex> globalFaces;
1167  faceList mergedFaces;
1168  pointField mergedPoints;
1170  (
1171  mesh,
1172  patch.localFaces(),
1173  patch.meshPoints(),
1174  patch.meshPointMap(),
1175  pointToGlobal,
1176  uniqueMeshPointLabels,
1177  globalPoints,
1178  globalFaces,
1179  mergedFaces,
1180  mergedPoints
1181  );
1182 
1183  // Collect the patch coverage
1184  scalarField mergedCoverage;
1185  globalFaces().gather
1186  (
1187  UPstream::worldComm,
1188  labelList(UPstream::procID(UPstream::worldComm)),
1189  patchCoverage[patchi],
1190  mergedCoverage
1191  );
1192 
1193  // Write the surface
1194  if (Pstream::master())
1195  {
1197  (
1198  mesh.time().writeFormat(),
1200  ).write
1201  (
1202  outputPath,
1203  patch.name() + "_coverage",
1204  mergedPoints,
1205  mergedFaces,
1206  false,
1207  "coverage",
1208  mergedCoverage
1209  );
1210  }
1211  }
1212  }
1213  }
1214 
1215  return noFailedChecks;
1216 }
1217 
1218 
1219 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
void resize(const label newSize)
Resize the hash table for efficiency.
Definition: HashTable.C:506
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:352
const word & name() const
Return name.
Definition: IOobject.H:307
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &p, Field< typename PrimitivePatch< FaceList, PointField >::PointType > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::FaceType > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
const Map< label > & meshPointMap() const
Mesh point map. Given the global point index find its.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const Field< PointType > & points() const
Return reference to global points.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
fileName globalPath() const
Return the global path.
Definition: TimePaths.H:132
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:281
IOstream::compressionType writeCompression() const
Default write compression.
Definition: Time.H:293
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:60
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:66
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:90
A collection of cell labels.
Definition: cellSet.H:51
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
virtual bool owner() const =0
Does this side own the patch ?
static scalarField calcFaceTol(const UList< face > &faces, const pointField &points, const pointField &faceCentres)
Calculate typical tolerance per face. Is currently max distance.
const word & name() const
Return const reference to name.
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:67
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
Definition: ops.H:70
A list of face labels.
Definition: faceSet.H:51
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
A class for handling file names.
Definition: fileName.H:82
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:101
Foam::intersection.
Definition: intersection.H:50
Non-conformal cyclic poly patch. As nonConformalCoupledPolyPatch, but the neighbouring patch is local...
virtual bool owner() const
Inherit the cyclic owner method.
const patchToPatches::intersection & intersection() const
Access the intersection engine.
const nonConformalCyclicPolyPatch & nbrPatch() const
Neighbour patch.
const polyPatch & origPatch() const
Original patch.
label index() const
Return the index of this patch in the boundaryMesh.
const word & name() const
Return name.
Class to generate patchToPatch coupling geometry. A full geometric intersection is done between a fac...
A set of point labels.
Definition: pointSet.H:51
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1023
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1040
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:988
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:407
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:1029
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1012
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:276
label nInternalFaces() const
label nCells() const
label nPoints() const
label nFaces() const
virtual bool write(const bool write=true) const
Write using setting from DB.
Base class for writing coordinate sets with data.
Definition: setWriter.H:64
Default transformation behaviour for position.
A surfaceWriter for VTK legacy format with support for writing ASCII or binary.
virtual void write(const fileName &outputDir, const fileName &surfaceName, const pointField &points, const faceList &faces, const wordList &fieldNames, const bool writePointValues #define FieldTypeValuesConstArg(Type, nullArg)) const
Write fields for a single surface to file.
Wedge front and back plane patch.
const vector & centreNormal() const
Return plane normal between the wedge boundaries.
scalar cosAngle() const
Return the cosine of the wedge angle.
const vector & n() const
Return the normal to the patch.
const vector & axis() const
Return axis of the wedge.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const pointField & points
label nPoints
const cellShapeList & cells
const fvPatchList & patches
Tools for checking the mesh.
Functions for checking mesh topology and geometry.
const dimensionedScalar mp
Proton mass.
const dimensionedScalar e
Elementary charge.
bool checkFaceWeight(const polyMesh &mesh, const bool report, const scalar minWeight=0.05, labelHashSet *setPtr=nullptr)
Check for face weights.
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.
bool checkEdgeAlignment(const polyMesh &mesh, const bool report, const Vector< label > &directions, labelHashSet *setPtr)
Check edge alignment for 1D/2D cases.
bool checkEdgeLength(const primitiveMesh &mesh, const bool report, const scalar minLenSqr, labelHashSet *setPtr=nullptr)
Check edge length.
bool checkCoupledPoints(const polyMesh &, const bool report, labelHashSet *)
Check 0th vertex on coupled faces.
bool checkConcaveCells(const primitiveMesh &mesh, const scalar planarCosAngle, const bool report=false, labelHashSet *setPtr=nullptr)
Check for concave cells by the planes of faces.
bool checkWedges(const polyMesh &, const bool report, const Vector< label > &, labelHashSet *)
Check wedge orientation.
Definition: checkGeometry.C:87
bool checkPointNearness(const primitiveMesh &mesh, const bool report, const scalar reportDistSqr, labelHashSet *setPtr=nullptr)
Check for point-point-nearness,.
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.
label findOppositeWedge(const polyMesh &, const wedgePolyPatch &)
Find wedge with opposite orientation. Note: does not actually check.
Definition: checkGeometry.C:47
bool checkClosedCells(const primitiveMesh &mesh, const scalar closedThreshold, const scalar aspectThreshold, const bool report=false, labelHashSet *setPtr=nullptr, labelHashSet *highAspectSetPtr=nullptr, const Vector< label > &solutionD=Vector< label >::one)
Check cells for closedness.
void mergeAndWrite(const polyMesh &mesh, const surfaceWriter &setWriter, const word &name, const indirectPrimitivePatch setPatch, const fileName &outputDir)
Generate merged surface on master and write. Needs input patch.
bool checkClosedBoundary(const primitiveMesh &mesh, const scalar closedThreshold, const bool report=false)
Check boundary for closedness.
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 checkCellVolumes(const primitiveMesh &mesh, const bool report=false, labelHashSet *setPtr=nullptr)
Check for negative cell volumes.
label checkGeometry(const polyMesh &mesh, const bool allGeometry, const scalar nonOrthThreshold, const scalar skewThreshold, const autoPtr< surfaceWriter > &, const autoPtr< setWriter > &)
Check the geometry.
bool checkFaceAreas(const primitiveMesh &mesh, const bool report=false, labelHashSet *setPtr=nullptr)
Check for negative face areas.
bool checkFaceSkewness(const polyMesh &mesh, const scalar skewThreshold, const bool report=false, labelHashSet *setPtr=nullptr)
Check face skewness.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
scalar radToDeg(const scalar rad)
Convert radians to degrees.
scalar degToRad(const scalar deg)
Convert degrees to radians.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:258
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
uint8_t direction
Definition: direction.H:45
dimensionedScalar acos(const dimensionedScalar &ds)
labelList f(nPoints)
volScalarField & p