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