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