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-2022 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 "PatchTools.H"
27 #include "checkGeometry.H"
28 #include "polyMesh.H"
29 #include "cellSet.H"
30 #include "faceSet.H"
31 #include "pointSet.H"
32 #include "EdgeMap.H"
33 #include "wedgePolyPatch.H"
34 #include "unitConversion.H"
36 
37 #include "vtkSurfaceWriter.H"
38 #include "setWriter.H"
39 #include "writeFile.H"
41 
42 #include "checkTools.H"
43 #include "Time.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
48 (
49  const polyMesh& mesh,
50  const wedgePolyPatch& wpp
51 )
52 {
53  const polyBoundaryMesh& patches = mesh.boundaryMesh();
54 
55  scalar wppCosAngle = wpp.cosAngle();
56 
58  {
59  if
60  (
61  patchi != wpp.index()
62  && patches[patchi].size()
63  && isA<wedgePolyPatch>(patches[patchi])
64  )
65  {
66  const wedgePolyPatch& pp =
67  refCast<const wedgePolyPatch>(patches[patchi]);
68 
69  // Calculate (cos of) angle to wpp (not pp!) centre normal
70  scalar ppCosAngle = wpp.centreNormal() & pp.n();
71 
72  if
73  (
74  pp.size() == wpp.size()
75  && mag(pp.axis() & wpp.axis()) >= (1-1e-3)
76  && mag(ppCosAngle - wppCosAngle) >= 1e-3
77  )
78  {
79  return patchi;
80  }
81  }
82  }
83  return -1;
84 }
85 
86 
88 (
89  const polyMesh& mesh,
90  const bool report,
91  const Vector<label>& directions,
92  labelHashSet* setPtr
93 )
94 {
95  // To mark edges without calculating edge addressing
96  EdgeMap<label> edgesInError;
97 
98  const pointField& p = mesh.points();
99  const faceList& fcs = mesh.faces();
100 
101 
102  const polyBoundaryMesh& patches = mesh.boundaryMesh();
104  {
105  if (patches[patchi].size() && isA<wedgePolyPatch>(patches[patchi]))
106  {
107  const wedgePolyPatch& pp =
108  refCast<const wedgePolyPatch>(patches[patchi]);
109 
110  scalar wedgeAngle = acos(pp.cosAngle());
111 
112  if (report)
113  {
114  Info<< " Wedge " << pp.name() << " with angle "
115  << radToDeg(wedgeAngle) << " degrees"
116  << endl;
117  }
118 
119  // Find opposite
120  label oppositePatchi = findOppositeWedge(mesh, pp);
121 
122  if (oppositePatchi == -1)
123  {
124  if (report)
125  {
126  Info<< " ***Cannot find opposite wedge for wedge "
127  << pp.name() << endl;
128  }
129  return true;
130  }
131 
132  const wedgePolyPatch& opp =
133  refCast<const wedgePolyPatch>(patches[oppositePatchi]);
134 
135 
136  if (mag(opp.axis() & pp.axis()) < (1-1e-3))
137  {
138  if (report)
139  {
140  Info<< " ***Wedges do not have the same axis."
141  << " Encountered " << pp.axis()
142  << " on patch " << pp.name()
143  << " which differs from " << opp.axis()
144  << " on opposite wedge patch" << opp.axis()
145  << endl;
146  }
147  return true;
148  }
149 
150 
151 
152  // Mark edges on wedgePatches
153  forAll(pp, i)
154  {
155  const face& f = pp[i];
156  forAll(f, fp)
157  {
158  label p0 = f[fp];
159  label p1 = f.nextLabel(fp);
160  edgesInError.insert(edge(p0, p1), -1); // non-error value
161  }
162  }
163 
164 
165  // Check that wedge patch is flat
166  const point& p0 = p[pp.meshPoints()[0]];
167  forAll(pp.meshPoints(), i)
168  {
169  const point& pt = p[pp.meshPoints()[i]];
170  scalar d = mag((pt - p0) & pp.n());
171 
172  if (d > rootSmall)
173  {
174  if (report)
175  {
176  Info<< " ***Wedge patch " << pp.name() << " not planar."
177  << " Point " << pt << " is not in patch plane by "
178  << d << " metre."
179  << endl;
180  }
181  return true;
182  }
183  }
184  }
185  }
186 
187 
188 
189  // Check all non-wedge faces
190  label nEdgesInError = 0;
191 
192  forAll(fcs, facei)
193  {
194  const face& f = fcs[facei];
195 
196  forAll(f, fp)
197  {
198  label p0 = f[fp];
199  label p1 = f.nextLabel(fp);
200  if (p0 < p1)
201  {
202  vector d(p[p1]-p[p0]);
203  scalar magD = mag(d);
204 
205  if (magD > rootVSmall)
206  {
207  d /= magD;
208 
209  // Check how many empty directions are used by the edge.
210  label nEmptyDirs = 0;
211  label nNonEmptyDirs = 0;
212  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
213  {
214  if (mag(d[cmpt]) > 1e-6)
215  {
216  if (directions[cmpt] == 0)
217  {
218  nEmptyDirs++;
219  }
220  else
221  {
222  nNonEmptyDirs++;
223  }
224  }
225  }
226 
227  if (nEmptyDirs == 0)
228  {
229  // Purely in ok directions.
230  }
231  else if (nEmptyDirs == 1)
232  {
233  // Ok if purely in empty directions.
234  if (nNonEmptyDirs > 0)
235  {
236  if (edgesInError.insert(edge(p0, p1), facei))
237  {
238  nEdgesInError++;
239  }
240  }
241  }
242  else if (nEmptyDirs > 1)
243  {
244  // Always an error
245  if (edgesInError.insert(edge(p0, p1), facei))
246  {
247  nEdgesInError++;
248  }
249  }
250  }
251  }
252  }
253  }
254 
255  label nErrorEdges = returnReduce(nEdgesInError, sumOp<label>());
256 
257  if (nErrorEdges > 0)
258  {
259  if (report)
260  {
261  Info<< " ***Number of edges not aligned with or perpendicular to "
262  << "non-empty directions: " << nErrorEdges << endl;
263  }
264 
265  if (setPtr)
266  {
267  setPtr->resize(2*nEdgesInError);
268  forAllConstIter(EdgeMap<label>, edgesInError, iter)
269  {
270  if (iter() >= 0)
271  {
272  setPtr->insert(iter.key()[0]);
273  setPtr->insert(iter.key()[1]);
274  }
275  }
276  }
277 
278  return true;
279  }
280  else
281  {
282  if (report)
283  {
284  Info<< " All edges aligned with or perpendicular to "
285  << "non-empty directions." << endl;
286  }
287  return false;
288  }
289 }
290 
291 
292 namespace Foam
293 {
294  //- Default transformation behaviour for position
295  class transformPositionList
296  {
297  public:
298 
299  //- Transform patch-based field
300  void operator()
301  (
302  const coupledPolyPatch& cpp,
303  List<pointField>& pts
304  ) const
305  {
306  // Each element of pts is all the points in the face. Convert into
307  // lists of size cpp to transform.
308 
309  List<pointField> newPts(pts.size());
310  forAll(pts, facei)
311  {
312  newPts[facei].setSize(pts[facei].size());
313  }
314 
315  label index = 0;
316  while (true)
317  {
318  label n = 0;
319 
320  // Extract for every face the i'th position
321  pointField ptsAtIndex(pts.size(), Zero);
322  forAll(cpp, facei)
323  {
324  const pointField& facePts = pts[facei];
325  if (facePts.size() > index)
326  {
327  ptsAtIndex[facei] = facePts[index];
328  n++;
329  }
330  }
331 
332  if (n == 0)
333  {
334  break;
335  }
336 
337  // Now ptsAtIndex will have for every face either zero or
338  // the position of the i'th vertex. Transform.
339  cpp.transform().transformPosition(ptsAtIndex, ptsAtIndex);
340 
341  // Extract back from ptsAtIndex into newPts
342  forAll(cpp, facei)
343  {
344  pointField& facePts = newPts[facei];
345  if (facePts.size() > index)
346  {
347  facePts[index] = ptsAtIndex[facei];
348  }
349  }
350 
351  index++;
352  }
353 
354  pts.transfer(newPts);
355  }
356  };
357 }
358 
359 
361 (
362  const polyMesh& mesh,
363  const bool report,
364  labelHashSet* setPtr
365 )
366 {
367  const pointField& p = mesh.points();
368  const faceList& fcs = mesh.faces();
369  const polyBoundaryMesh& patches = mesh.boundaryMesh();
370 
371  // Zero'th point on coupled faces
372  // pointField nbrZeroPoint(fcs.size()-mesh.nInternalFaces(), vector::max);
373  List<pointField> nbrPoints(fcs.size() - mesh.nInternalFaces());
374 
375  // Exchange zero point
377  {
378  if (patches[patchi].coupled())
379  {
380  const coupledPolyPatch& cpp = refCast<const coupledPolyPatch>
381  (
382  patches[patchi]
383  );
384 
385  forAll(cpp, i)
386  {
387  label bFacei = cpp.start() + i - mesh.nInternalFaces();
388  const face& f = cpp[i];
389  nbrPoints[bFacei].setSize(f.size());
390  forAll(f, fp)
391  {
392  const point& p0 = p[f[fp]];
393  nbrPoints[bFacei][fp] = p0;
394  }
395  }
396  }
397  }
398  syncTools::syncBoundaryFaceList
399  (
400  mesh,
401  nbrPoints,
402  eqOp<pointField>(),
403  transformPositionList()
404  );
405 
406  // Compare to local ones. Use same tolerance as for matching
407  label nErrorFaces = 0;
408  scalar avgMismatch = 0;
409  label nCoupledPoints = 0;
410 
412  {
413  if (patches[patchi].coupled())
414  {
415  const coupledPolyPatch& cpp =
416  refCast<const coupledPolyPatch>(patches[patchi]);
417 
418  if (cpp.owner())
419  {
420  scalarField smallDist
421  (
422  cpp.calcFaceTol
423  (
424  // cpp.matchTolerance(),
425  cpp,
426  cpp.points(),
427  cpp.faceCentres()
428  )
429  );
430 
431  forAll(cpp, i)
432  {
433  label bFacei = cpp.start() + i - mesh.nInternalFaces();
434  const face& f = cpp[i];
435 
436  if (f.size() != nbrPoints[bFacei].size())
437  {
439  << "Local face size : " << f.size()
440  << " does not equal neighbour face size : "
441  << nbrPoints[bFacei].size()
442  << abort(FatalError);
443  }
444 
445  label fp = 0;
446  forAll(f, j)
447  {
448  const point& p0 = p[f[fp]];
449  scalar d = mag(p0 - nbrPoints[bFacei][j]);
450 
451  if (d > smallDist[i])
452  {
453  if (setPtr)
454  {
455  setPtr->insert(cpp.start()+i);
456  }
457  nErrorFaces++;
458 
459  break;
460  }
461 
462  avgMismatch += d;
463  nCoupledPoints++;
464 
465  fp = f.rcIndex(fp);
466  }
467  }
468  }
469  }
470  }
471 
472  reduce(nErrorFaces, sumOp<label>());
473  reduce(avgMismatch, maxOp<scalar>());
474  reduce(nCoupledPoints, sumOp<label>());
475 
476  if (nCoupledPoints > 0)
477  {
478  avgMismatch /= nCoupledPoints;
479  }
480 
481  if (nErrorFaces > 0)
482  {
483  if (report)
484  {
485  Info<< " **Error in coupled point location: "
486  << nErrorFaces
487  << " faces have their 0th or consecutive vertex not opposite"
488  << " their coupled equivalent. Average mismatch "
489  << avgMismatch << "."
490  << endl;
491  }
492 
493  return true;
494  }
495  else
496  {
497  if (report)
498  {
499  Info<< " Coupled point location match (average "
500  << avgMismatch << ") OK." << endl;
501  }
502 
503  return false;
504  }
505 }
506 
507 
509 (
510  const polyMesh& mesh,
511  const bool allGeometry,
512  const autoPtr<surfaceWriter>& surfWriter,
513  const autoPtr<Foam::setWriter>& setWriter
514 )
515 {
516  label noFailedChecks = 0;
517 
518  Info<< "\nChecking geometry..." << endl;
519 
520  // Get a small relative length from the bounding box
521  const boundBox& globalBb = mesh.bounds();
522 
523  Info<< " Overall domain bounding box "
524  << globalBb.min() << " " << globalBb.max() << endl;
525 
526 
527  // Min length
528  scalar minDistSqr = magSqr(1e-6 * globalBb.span());
529 
530  // Geometric directions
531  const Vector<label> validDirs = (mesh.geometricD() + Vector<label>::one)/2;
532  Info<< " Mesh has " << mesh.nGeometricD()
533  << " geometric (non-empty/wedge) directions " << validDirs << endl;
534 
535  // Solution directions
536  const Vector<label> solDirs = (mesh.solutionD() + Vector<label>::one)/2;
537  Info<< " Mesh has " << mesh.nSolutionD()
538  << " solution (non-empty) directions " << solDirs << endl;
539 
540  if (mesh.nGeometricD() < 3)
541  {
542  pointSet nonAlignedPoints(mesh, "nonAlignedEdges", mesh.nPoints()/100);
543 
544  if
545  (
546  (
547  validDirs != solDirs
548  && checkWedges(mesh, true, validDirs, &nonAlignedPoints)
549  )
550  || (
551  validDirs == solDirs
552  && mesh.checkEdgeAlignment(true, validDirs, &nonAlignedPoints)
553  )
554  )
555  {
556  noFailedChecks++;
557  label nNonAligned = returnReduce
558  (
559  nonAlignedPoints.size(),
560  sumOp<label>()
561  );
562 
563  if (nNonAligned > 0)
564  {
565  Info<< " <<Writing " << nNonAligned
566  << " points on non-aligned edges to set "
567  << nonAlignedPoints.name() << endl;
568  nonAlignedPoints.instance() = mesh.pointsInstance();
569  nonAlignedPoints.write();
570  if (setWriter.valid())
571  {
572  mergeAndWrite(setWriter, nonAlignedPoints);
573  }
574  }
575  }
576  }
577 
578  if (mesh.checkClosedBoundary(true)) noFailedChecks++;
579 
580  {
581  cellSet cells(mesh, "nonClosedCells", mesh.nCells()/100+1);
582  cellSet aspectCells(mesh, "highAspectRatioCells", mesh.nCells()/100+1);
583  if
584  (
585  mesh.checkClosedCells
586  (
587  true,
588  &cells,
589  &aspectCells,
590  mesh.geometricD()
591  )
592  )
593  {
594  noFailedChecks++;
595 
596  label nNonClosed = returnReduce(cells.size(), sumOp<label>());
597 
598  if (nNonClosed > 0)
599  {
600  Info<< " <<Writing " << nNonClosed
601  << " non closed cells to set " << cells.name() << endl;
602  cells.instance() = mesh.pointsInstance();
603  cells.write();
604  if (surfWriter.valid())
605  {
606  mergeAndWrite(surfWriter(), cells);
607  }
608  }
609  }
610 
611  label nHighAspect = returnReduce(aspectCells.size(), sumOp<label>());
612 
613  if (nHighAspect > 0)
614  {
615  Info<< " <<Writing " << nHighAspect
616  << " cells with high aspect ratio to set "
617  << aspectCells.name() << endl;
618  aspectCells.instance() = mesh.pointsInstance();
619  aspectCells.write();
620  if (surfWriter.valid())
621  {
622  mergeAndWrite(surfWriter(), aspectCells);
623  }
624  }
625  }
626 
627  {
628  faceSet faces(mesh, "zeroAreaFaces", mesh.nFaces()/100+1);
629  if (mesh.checkFaceAreas(true, &faces))
630  {
631  noFailedChecks++;
632 
633  label nFaces = returnReduce(faces.size(), sumOp<label>());
634 
635  if (nFaces > 0)
636  {
637  Info<< " <<Writing " << nFaces
638  << " zero area faces to set " << faces.name() << endl;
639  faces.instance() = mesh.pointsInstance();
640  faces.write();
641  if (surfWriter.valid())
642  {
643  mergeAndWrite(surfWriter(), faces);
644  }
645  }
646  }
647  }
648 
649  {
650  cellSet cells(mesh, "zeroVolumeCells", mesh.nCells()/100+1);
651  if (mesh.checkCellVolumes(true, &cells))
652  {
653  noFailedChecks++;
654 
655  label nCells = returnReduce(cells.size(), sumOp<label>());
656 
657  if (nCells > 0)
658  {
659  Info<< " <<Writing " << nCells
660  << " zero volume cells to set " << cells.name() << endl;
661  cells.instance() = mesh.pointsInstance();
662  cells.write();
663  if (surfWriter.valid())
664  {
665  mergeAndWrite(surfWriter(), cells);
666  }
667  }
668  }
669  }
670 
671  {
672  faceSet faces(mesh, "nonOrthoFaces", mesh.nFaces()/100+1);
673  if (mesh.checkFaceOrthogonality(true, &faces))
674  {
675  noFailedChecks++;
676  }
677 
678  label nFaces = returnReduce(faces.size(), sumOp<label>());
679 
680  if (nFaces > 0)
681  {
682  Info<< " <<Writing " << nFaces
683  << " non-orthogonal faces to set " << faces.name() << endl;
684  faces.instance() = mesh.pointsInstance();
685  faces.write();
686  if (surfWriter.valid())
687  {
688  mergeAndWrite(surfWriter(), faces);
689  }
690  }
691  }
692 
693  {
694  faceSet faces(mesh, "wrongOrientedFaces", mesh.nFaces()/100 + 1);
695  if (mesh.checkFacePyramids(true, -small, &faces))
696  {
697  noFailedChecks++;
698 
699  label nFaces = returnReduce(faces.size(), sumOp<label>());
700 
701  if (nFaces > 0)
702  {
703  Info<< " <<Writing " << nFaces
704  << " faces with incorrect orientation to set "
705  << faces.name() << endl;
706  faces.instance() = mesh.pointsInstance();
707  faces.write();
708  if (surfWriter.valid())
709  {
710  mergeAndWrite(surfWriter(), faces);
711  }
712  }
713  }
714  }
715 
716  {
717  faceSet faces(mesh, "skewFaces", mesh.nFaces()/100+1);
718  if (mesh.checkFaceSkewness(true, &faces))
719  {
720  noFailedChecks++;
721 
722  label nFaces = returnReduce(faces.size(), sumOp<label>());
723 
724  if (nFaces > 0)
725  {
726  Info<< " <<Writing " << nFaces
727  << " skew faces to set " << faces.name() << endl;
728  faces.instance() = mesh.pointsInstance();
729  faces.write();
730  if (surfWriter.valid())
731  {
732  mergeAndWrite(surfWriter(), faces);
733  }
734  }
735  }
736  }
737 
738  {
739  faceSet faces(mesh, "coupledFaces", mesh.nFaces()/100 + 1);
740  if (checkCoupledPoints(mesh, true, &faces))
741  {
742  noFailedChecks++;
743 
744  label nFaces = returnReduce(faces.size(), sumOp<label>());
745 
746  if (nFaces > 0)
747  {
748  Info<< " <<Writing " << nFaces
749  << " faces with incorrectly matched 0th (or consecutive)"
750  << " vertex to set "
751  << faces.name() << endl;
752  faces.instance() = mesh.pointsInstance();
753  faces.write();
754  if (surfWriter.valid())
755  {
756  mergeAndWrite(surfWriter(), faces);
757  }
758  }
759  }
760  }
761 
762  if (allGeometry)
763  {
764  faceSet faces(mesh, "lowQualityTetFaces", mesh.nFaces()/100+1);
765  if
766  (
768  (
769  mesh,
770  polyMeshTetDecomposition::minTetQuality,
771  true,
772  &faces
773  )
774  )
775  {
776  noFailedChecks++;
777 
778  label nFaces = returnReduce(faces.size(), sumOp<label>());
779 
780  if (nFaces > 0)
781  {
782  Info<< " <<Writing " << nFaces
783  << " faces with low quality or negative volume "
784  << "decomposition tets to set " << faces.name() << endl;
785  faces.instance() = mesh.pointsInstance();
786  faces.write();
787  if (surfWriter.valid())
788  {
789  mergeAndWrite(surfWriter(), faces);
790  }
791  }
792  }
793  }
794 
795  if (allGeometry)
796  {
797  // Note use of nPoints since don't want edge construction.
798  pointSet points(mesh, "shortEdges", mesh.nPoints()/1000 + 1);
799  if (mesh.checkEdgeLength(true, minDistSqr, &points))
800  {
801  // noFailedChecks++;
802 
803  label nPoints = returnReduce(points.size(), sumOp<label>());
804 
805  if (nPoints > 0)
806  {
807  Info<< " <<Writing " << nPoints
808  << " points on short edges to set " << points.name()
809  << endl;
810  points.instance() = mesh.pointsInstance();
811  points.write();
812  if (setWriter.valid())
813  {
814  mergeAndWrite(setWriter, points);
815  }
816  }
817  }
818 
819  label nEdgeClose = returnReduce(points.size(), sumOp<label>());
820 
821  if (mesh.checkPointNearness(false, minDistSqr, &points))
822  {
823  // noFailedChecks++;
824 
825  label nPoints = returnReduce(points.size(), sumOp<label>());
826 
827  if (nPoints > nEdgeClose)
828  {
829  pointSet nearPoints(mesh, "nearPoints", points);
830  Info<< " <<Writing " << nPoints
831  << " near (closer than " << Foam::sqrt(minDistSqr)
832  << " apart) points to set " << nearPoints.name() << endl;
833  nearPoints.instance() = mesh.pointsInstance();
834  nearPoints.write();
835  if (setWriter.valid())
836  {
837  mergeAndWrite(setWriter, nearPoints);
838  }
839  }
840  }
841  }
842 
843  if (allGeometry)
844  {
845  faceSet faces(mesh, "concaveFaces", mesh.nFaces()/100 + 1);
846  if (mesh.checkFaceAngles(true, 10, &faces))
847  {
848  // noFailedChecks++;
849 
850  label nFaces = returnReduce(faces.size(), sumOp<label>());
851 
852  if (nFaces > 0)
853  {
854  Info<< " <<Writing " << nFaces
855  << " faces with concave angles to set " << faces.name()
856  << endl;
857  faces.instance() = mesh.pointsInstance();
858  faces.write();
859  if (surfWriter.valid())
860  {
861  mergeAndWrite(surfWriter(), faces);
862  }
863  }
864  }
865  }
866 
867  if (allGeometry)
868  {
869  faceSet faces(mesh, "warpedFaces", mesh.nFaces()/100 + 1);
870  if (mesh.checkFaceFlatness(true, 0.8, &faces))
871  {
872  // noFailedChecks++;
873 
874  label nFaces = returnReduce(faces.size(), sumOp<label>());
875 
876  if (nFaces > 0)
877  {
878  Info<< " <<Writing " << nFaces
879  << " warped faces to set " << faces.name() << endl;
880  faces.instance() = mesh.pointsInstance();
881  faces.write();
882  if (surfWriter.valid())
883  {
884  mergeAndWrite(surfWriter(), faces);
885  }
886  }
887  }
888  }
889 
890  if (allGeometry)
891  {
892  cellSet cells(mesh, "underdeterminedCells", mesh.nCells()/100);
893  if (mesh.checkCellDeterminant(true, &cells))
894  {
895  noFailedChecks++;
896 
897  label nCells = returnReduce(cells.size(), sumOp<label>());
898 
899  Info<< " <<Writing " << nCells
900  << " under-determined cells to set " << cells.name() << endl;
901  cells.instance() = mesh.pointsInstance();
902  cells.write();
903  if (surfWriter.valid())
904  {
905  mergeAndWrite(surfWriter(), cells);
906  }
907  }
908  }
909 
910  if (allGeometry)
911  {
912  cellSet cells(mesh, "concaveCells", mesh.nCells()/100);
913  if (mesh.checkConcaveCells(true, &cells))
914  {
915  noFailedChecks++;
916 
917  label nCells = returnReduce(cells.size(), sumOp<label>());
918 
919  Info<< " <<Writing " << nCells
920  << " concave cells to set " << cells.name() << endl;
921  cells.instance() = mesh.pointsInstance();
922  cells.write();
923  if (surfWriter.valid())
924  {
925  mergeAndWrite(surfWriter(), cells);
926  }
927  }
928  }
929 
930  if (allGeometry)
931  {
932  faceSet faces(mesh, "lowWeightFaces", mesh.nFaces()/100);
933  if (mesh.checkFaceWeight(true, 0.05, &faces))
934  {
935  noFailedChecks++;
936 
937  label nFaces = returnReduce(faces.size(), sumOp<label>());
938 
939  Info<< " <<Writing " << nFaces
940  << " faces with low interpolation weights to set "
941  << faces.name() << endl;
942  faces.instance() = mesh.pointsInstance();
943  faces.write();
944  if (surfWriter.valid())
945  {
946  mergeAndWrite(surfWriter(), faces);
947  }
948  }
949  }
950 
951  if (allGeometry)
952  {
953  faceSet faces(mesh, "lowVolRatioFaces", mesh.nFaces()/100);
954  if (mesh.checkVolRatio(true, 0.01, &faces))
955  {
956  noFailedChecks++;
957 
958  label nFaces = returnReduce(faces.size(), sumOp<label>());
959 
960  Info<< " <<Writing " << nFaces
961  << " faces with low volume ratio cells to set "
962  << faces.name() << endl;
963  faces.instance() = mesh.pointsInstance();
964  faces.write();
965  if (surfWriter.valid())
966  {
967  mergeAndWrite(surfWriter(), faces);
968  }
969  }
970  }
971 
972  if (allGeometry)
973  {
974  const fileName outputPath =
975  mesh.time().globalPath()
976  /functionObjects::writeFile::outputPrefix
977  /(mesh.name() != polyMesh::defaultRegion ? mesh.name() : word())
978  /"checkMesh"
979  /mesh.time().name();
980 
981  const polyBoundaryMesh& patches = mesh.boundaryMesh();
982 
983  // Compute coverage for all orig patches
984  PtrList<scalarField> patchCoverage(patches.size());
985  forAll(patches, nccPatchi)
986  {
987  if (isA<nonConformalCyclicPolyPatch>(patches[nccPatchi]))
988  {
989  const nonConformalCyclicPolyPatch& nccPp =
990  refCast<const nonConformalCyclicPolyPatch>
991  (patches[nccPatchi]);
992 
993  if (nccPp.owner())
994  {
995  const polyPatch& origPp = nccPp.origPatch();
996  const polyPatch& nbrOrigPp = nccPp.nbrPatch().origPatch();
997 
998  const patchToPatches::intersection& intersection =
999  nccPp.intersection();
1000 
1001  if (!patchCoverage.set(origPp.index()))
1002  {
1003  patchCoverage.set
1004  (
1005  origPp.index(),
1006  scalarField(origPp.size(), 0)
1007  );
1008  }
1009 
1010  patchCoverage[origPp.index()] +=
1011  intersection.srcCoverage();
1012 
1013  if (!patchCoverage.set(nbrOrigPp.index()))
1014  {
1015  patchCoverage.set
1016  (
1017  nbrOrigPp.index(),
1018  scalarField(nbrOrigPp.size(), 0)
1019  );
1020  }
1021 
1022  patchCoverage[nbrOrigPp.index()] +=
1023  intersection.tgtCoverage();
1024  }
1025  }
1026  }
1027 
1028  // Write out to surface files
1030  {
1031  if (patchCoverage.set(patchi))
1032  {
1033  const polyPatch& patch = patches[patchi];
1034 
1035  // Collect the patch geometry
1036  labelList pointToGlobal;
1037  labelList uniqueMeshPointLabels;
1038  autoPtr<globalIndex> globalPoints;
1039  autoPtr<globalIndex> globalFaces;
1040  faceList mergedFaces;
1041  pointField mergedPoints;
1043  (
1044  mesh,
1045  patch.localFaces(),
1046  patch.meshPoints(),
1047  patch.meshPointMap(),
1048  pointToGlobal,
1049  uniqueMeshPointLabels,
1050  globalPoints,
1051  globalFaces,
1052  mergedFaces,
1053  mergedPoints
1054  );
1055 
1056  // Collect the patch coverage
1057  scalarField mergedCoverage;
1058  globalFaces().gather
1059  (
1060  UPstream::worldComm,
1061  labelList(UPstream::procID(UPstream::worldComm)),
1062  patchCoverage[patchi],
1063  mergedCoverage
1064  );
1065 
1066  // Write the surface
1067  if (Pstream::master())
1068  {
1069  vtkSurfaceWriter
1070  (
1071  mesh.time().writeFormat(),
1072  mesh.time().writeCompression()
1073  ).write
1074  (
1075  outputPath,
1076  patch.name() + "_coverage",
1077  mergedPoints,
1078  mergedFaces,
1079  false,
1080  "coverage",
1081  mergedCoverage
1082  );
1083  }
1084  }
1085  }
1086  }
1087 
1088  return noFailedChecks;
1089 }
1090 
1091 
1092 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Routines for checking mesh geometry.
Tools for checking the mesh.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
virtual Ostream & write(const char)=0
Write character.
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.
static const Form one
Definition: VectorSpace.H:114
const word & name() const
Return const reference to name.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
volScalarField scalarField(fieldObject, mesh)
label patchi
const pointField & points
label nPoints
const cellShapeList & cells
const fvPatchList & patches
const dimensionedScalar e
Elementary charge.
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.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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
label findOppositeWedge(const polyMesh &, const wedgePolyPatch &)
Find wedge with opposite orientation. Note: does not actually check.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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.
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
vector point
Point is a vector.
Definition: point.H:41
label checkGeometry(const polyMesh &mesh, const bool allGeometry, const autoPtr< surfaceWriter > &, const autoPtr< setWriter > &)
Check the geometry.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
bool checkWedges(const polyMesh &, const bool report, const Vector< label > &, labelHashSet *)
Check wedge orientation.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
error FatalError
bool checkCoupledPoints(const polyMesh &, const bool report, labelHashSet *)
Check 0th vertex on coupled faces.
List< face > faceList
Definition: faceListFwd.H:43
dimensioned< scalar > magSqr(const dimensioned< Type > &)
uint8_t direction
Definition: direction.H:45
dimensionedScalar acos(const dimensionedScalar &ds)
labelList f(nPoints)
volScalarField & p
Unit conversion functions.