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 "writer.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.transformPosition(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 bool allGeometry,
486  const autoPtr<surfaceWriter>& surfWriter,
487  const autoPtr<writer<scalar>>& setWriter
488 )
489 {
490  label noFailedChecks = 0;
491 
492  Info<< "\nChecking geometry..." << endl;
493 
494  // Get a small relative length from the bounding box
495  const boundBox& globalBb = mesh.bounds();
496 
497  Info<< " Overall domain bounding box "
498  << globalBb.min() << " " << globalBb.max() << endl;
499 
500 
501  // Min length
502  scalar minDistSqr = magSqr(1e-6 * globalBb.span());
503 
504  // Geometric directions
505  const Vector<label> validDirs = (mesh.geometricD() + Vector<label>::one)/2;
506  Info<< " Mesh has " << mesh.nGeometricD()
507  << " geometric (non-empty/wedge) directions " << validDirs << endl;
508 
509  // Solution directions
510  const Vector<label> solDirs = (mesh.solutionD() + Vector<label>::one)/2;
511  Info<< " Mesh has " << mesh.nSolutionD()
512  << " solution (non-empty) directions " << solDirs << endl;
513 
514  if (mesh.nGeometricD() < 3)
515  {
516  pointSet nonAlignedPoints(mesh, "nonAlignedEdges", mesh.nPoints()/100);
517 
518  if
519  (
520  (
521  validDirs != solDirs
522  && checkWedges(mesh, true, validDirs, &nonAlignedPoints)
523  )
524  || (
525  validDirs == solDirs
526  && mesh.checkEdgeAlignment(true, validDirs, &nonAlignedPoints)
527  )
528  )
529  {
530  noFailedChecks++;
531  label nNonAligned = returnReduce
532  (
533  nonAlignedPoints.size(),
534  sumOp<label>()
535  );
536 
537  if (nNonAligned > 0)
538  {
539  Info<< " <<Writing " << nNonAligned
540  << " points on non-aligned edges to set "
541  << nonAlignedPoints.name() << endl;
542  nonAlignedPoints.instance() = mesh.pointsInstance();
543  nonAlignedPoints.write();
544  if (setWriter.valid())
545  {
546  mergeAndWrite(setWriter, nonAlignedPoints);
547  }
548  }
549  }
550  }
551 
552  if (mesh.checkClosedBoundary(true)) noFailedChecks++;
553 
554  {
555  cellSet cells(mesh, "nonClosedCells", mesh.nCells()/100+1);
556  cellSet aspectCells(mesh, "highAspectRatioCells", mesh.nCells()/100+1);
557  if
558  (
559  mesh.checkClosedCells
560  (
561  true,
562  &cells,
563  &aspectCells,
564  mesh.geometricD()
565  )
566  )
567  {
568  noFailedChecks++;
569 
570  label nNonClosed = returnReduce(cells.size(), sumOp<label>());
571 
572  if (nNonClosed > 0)
573  {
574  Info<< " <<Writing " << nNonClosed
575  << " non closed cells to set " << cells.name() << endl;
576  cells.instance() = mesh.pointsInstance();
577  cells.write();
578  if (surfWriter.valid())
579  {
580  mergeAndWrite(surfWriter(), cells);
581  }
582  }
583  }
584 
585  label nHighAspect = returnReduce(aspectCells.size(), sumOp<label>());
586 
587  if (nHighAspect > 0)
588  {
589  Info<< " <<Writing " << nHighAspect
590  << " cells with high aspect ratio to set "
591  << aspectCells.name() << endl;
592  aspectCells.instance() = mesh.pointsInstance();
593  aspectCells.write();
594  if (surfWriter.valid())
595  {
596  mergeAndWrite(surfWriter(), aspectCells);
597  }
598  }
599  }
600 
601  {
602  faceSet faces(mesh, "zeroAreaFaces", mesh.nFaces()/100+1);
603  if (mesh.checkFaceAreas(true, &faces))
604  {
605  noFailedChecks++;
606 
607  label nFaces = returnReduce(faces.size(), sumOp<label>());
608 
609  if (nFaces > 0)
610  {
611  Info<< " <<Writing " << nFaces
612  << " zero area faces to set " << faces.name() << endl;
613  faces.instance() = mesh.pointsInstance();
614  faces.write();
615  if (surfWriter.valid())
616  {
617  mergeAndWrite(surfWriter(), faces);
618  }
619  }
620  }
621  }
622 
623  {
624  cellSet cells(mesh, "zeroVolumeCells", mesh.nCells()/100+1);
625  if (mesh.checkCellVolumes(true, &cells))
626  {
627  noFailedChecks++;
628 
629  label nCells = returnReduce(cells.size(), sumOp<label>());
630 
631  if (nCells > 0)
632  {
633  Info<< " <<Writing " << nCells
634  << " zero volume cells to set " << cells.name() << endl;
635  cells.instance() = mesh.pointsInstance();
636  cells.write();
637  if (surfWriter.valid())
638  {
639  mergeAndWrite(surfWriter(), cells);
640  }
641  }
642  }
643  }
644 
645  {
646  faceSet faces(mesh, "nonOrthoFaces", mesh.nFaces()/100+1);
647  if (mesh.checkFaceOrthogonality(true, &faces))
648  {
649  noFailedChecks++;
650  }
651 
652  label nFaces = returnReduce(faces.size(), sumOp<label>());
653 
654  if (nFaces > 0)
655  {
656  Info<< " <<Writing " << nFaces
657  << " non-orthogonal faces to set " << faces.name() << endl;
658  faces.instance() = mesh.pointsInstance();
659  faces.write();
660  if (surfWriter.valid())
661  {
662  mergeAndWrite(surfWriter(), faces);
663  }
664  }
665  }
666 
667  {
668  faceSet faces(mesh, "wrongOrientedFaces", mesh.nFaces()/100 + 1);
669  if (mesh.checkFacePyramids(true, -SMALL, &faces))
670  {
671  noFailedChecks++;
672 
673  label nFaces = returnReduce(faces.size(), sumOp<label>());
674 
675  if (nFaces > 0)
676  {
677  Info<< " <<Writing " << nFaces
678  << " faces with incorrect orientation to set "
679  << faces.name() << endl;
680  faces.instance() = mesh.pointsInstance();
681  faces.write();
682  if (surfWriter.valid())
683  {
684  mergeAndWrite(surfWriter(), faces);
685  }
686  }
687  }
688  }
689 
690  {
691  faceSet faces(mesh, "skewFaces", mesh.nFaces()/100+1);
692  if (mesh.checkFaceSkewness(true, &faces))
693  {
694  noFailedChecks++;
695 
696  label nFaces = returnReduce(faces.size(), sumOp<label>());
697 
698  if (nFaces > 0)
699  {
700  Info<< " <<Writing " << nFaces
701  << " skew faces to set " << faces.name() << endl;
702  faces.instance() = mesh.pointsInstance();
703  faces.write();
704  if (surfWriter.valid())
705  {
706  mergeAndWrite(surfWriter(), faces);
707  }
708  }
709  }
710  }
711 
712  {
713  faceSet faces(mesh, "coupledFaces", mesh.nFaces()/100 + 1);
714  if (checkCoupledPoints(mesh, true, &faces))
715  {
716  noFailedChecks++;
717 
718  label nFaces = returnReduce(faces.size(), sumOp<label>());
719 
720  if (nFaces > 0)
721  {
722  Info<< " <<Writing " << nFaces
723  << " faces with incorrectly matched 0th (or consecutive)"
724  << " vertex to set "
725  << faces.name() << endl;
726  faces.instance() = mesh.pointsInstance();
727  faces.write();
728  if (surfWriter.valid())
729  {
730  mergeAndWrite(surfWriter(), faces);
731  }
732  }
733  }
734  }
735 
736  if (allGeometry)
737  {
738  faceSet faces(mesh, "lowQualityTetFaces", mesh.nFaces()/100+1);
739  if
740  (
742  (
743  mesh,
745  true,
746  &faces
747  )
748  )
749  {
750  noFailedChecks++;
751 
752  label nFaces = returnReduce(faces.size(), sumOp<label>());
753 
754  if (nFaces > 0)
755  {
756  Info<< " <<Writing " << nFaces
757  << " faces with low quality or negative volume "
758  << "decomposition tets 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  if (allGeometry)
770  {
771  // Note use of nPoints since don't want edge construction.
772  pointSet points(mesh, "shortEdges", mesh.nPoints()/1000 + 1);
773  if (mesh.checkEdgeLength(true, minDistSqr, &points))
774  {
775  //noFailedChecks++;
776 
777  label nPoints = returnReduce(points.size(), sumOp<label>());
778 
779  if (nPoints > 0)
780  {
781  Info<< " <<Writing " << nPoints
782  << " points on short edges to set " << points.name()
783  << endl;
784  points.instance() = mesh.pointsInstance();
785  points.write();
786  if (setWriter.valid())
787  {
788  mergeAndWrite(setWriter, points);
789  }
790  }
791  }
792 
793  label nEdgeClose = returnReduce(points.size(), sumOp<label>());
794 
795  if (mesh.checkPointNearness(false, minDistSqr, &points))
796  {
797  //noFailedChecks++;
798 
799  label nPoints = returnReduce(points.size(), sumOp<label>());
800 
801  if (nPoints > nEdgeClose)
802  {
803  pointSet nearPoints(mesh, "nearPoints", points);
804  Info<< " <<Writing " << nPoints
805  << " near (closer than " << Foam::sqrt(minDistSqr)
806  << " apart) points to set " << nearPoints.name() << endl;
807  nearPoints.instance() = mesh.pointsInstance();
808  nearPoints.write();
809  if (setWriter.valid())
810  {
811  mergeAndWrite(setWriter, nearPoints);
812  }
813  }
814  }
815  }
816 
817  if (allGeometry)
818  {
819  faceSet faces(mesh, "concaveFaces", mesh.nFaces()/100 + 1);
820  if (mesh.checkFaceAngles(true, 10, &faces))
821  {
822  //noFailedChecks++;
823 
824  label nFaces = returnReduce(faces.size(), sumOp<label>());
825 
826  if (nFaces > 0)
827  {
828  Info<< " <<Writing " << nFaces
829  << " faces with concave angles to set " << faces.name()
830  << endl;
831  faces.instance() = mesh.pointsInstance();
832  faces.write();
833  if (surfWriter.valid())
834  {
835  mergeAndWrite(surfWriter(), faces);
836  }
837  }
838  }
839  }
840 
841  if (allGeometry)
842  {
843  faceSet faces(mesh, "warpedFaces", mesh.nFaces()/100 + 1);
844  if (mesh.checkFaceFlatness(true, 0.8, &faces))
845  {
846  //noFailedChecks++;
847 
848  label nFaces = returnReduce(faces.size(), sumOp<label>());
849 
850  if (nFaces > 0)
851  {
852  Info<< " <<Writing " << nFaces
853  << " warped faces to set " << faces.name() << endl;
854  faces.instance() = mesh.pointsInstance();
855  faces.write();
856  if (surfWriter.valid())
857  {
858  mergeAndWrite(surfWriter(), faces);
859  }
860  }
861  }
862  }
863 
864  if (allGeometry)
865  {
866  cellSet cells(mesh, "underdeterminedCells", mesh.nCells()/100);
867  if (mesh.checkCellDeterminant(true, &cells))
868  {
869  noFailedChecks++;
870 
871  label nCells = returnReduce(cells.size(), sumOp<label>());
872 
873  Info<< " <<Writing " << nCells
874  << " under-determined cells to set " << cells.name() << endl;
875  cells.instance() = mesh.pointsInstance();
876  cells.write();
877  if (surfWriter.valid())
878  {
879  mergeAndWrite(surfWriter(), cells);
880  }
881  }
882  }
883 
884  if (allGeometry)
885  {
886  cellSet cells(mesh, "concaveCells", mesh.nCells()/100);
887  if (mesh.checkConcaveCells(true, &cells))
888  {
889  noFailedChecks++;
890 
891  label nCells = returnReduce(cells.size(), sumOp<label>());
892 
893  Info<< " <<Writing " << nCells
894  << " concave cells to set " << cells.name() << endl;
895  cells.instance() = mesh.pointsInstance();
896  cells.write();
897  if (surfWriter.valid())
898  {
899  mergeAndWrite(surfWriter(), cells);
900  }
901  }
902  }
903 
904  if (allGeometry)
905  {
906  faceSet faces(mesh, "lowWeightFaces", mesh.nFaces()/100);
907  if (mesh.checkFaceWeight(true, 0.05, &faces))
908  {
909  noFailedChecks++;
910 
911  label nFaces = returnReduce(faces.size(), sumOp<label>());
912 
913  Info<< " <<Writing " << nFaces
914  << " faces with low interpolation weights to set "
915  << faces.name() << endl;
916  faces.instance() = mesh.pointsInstance();
917  faces.write();
918  if (surfWriter.valid())
919  {
920  mergeAndWrite(surfWriter(), faces);
921  }
922  }
923  }
924 
925  if (allGeometry)
926  {
927  faceSet faces(mesh, "lowVolRatioFaces", mesh.nFaces()/100);
928  if (mesh.checkVolRatio(true, 0.01, &faces))
929  {
930  noFailedChecks++;
931 
932  label nFaces = returnReduce(faces.size(), sumOp<label>());
933 
934  Info<< " <<Writing " << nFaces
935  << " faces with low volume ratio cells to set "
936  << faces.name() << endl;
937  faces.instance() = mesh.pointsInstance();
938  faces.write();
939  if (surfWriter.valid())
940  {
941  mergeAndWrite(surfWriter(), faces);
942  }
943  }
944  }
945 
946  if (allGeometry)
947  {
948  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
949 
950  const word tmName(mesh.time().timeName());
951  const word procAndTime(Foam::name(Pstream::myProcNo()) + "_" + tmName);
952 
953  autoPtr<surfaceWriter> patchWriter;
954  if (!surfWriter.valid())
955  {
956  patchWriter.reset(new vtkSurfaceWriter());
957  }
958  const surfaceWriter& wr =
959  (
960  surfWriter.valid()
961  ? surfWriter()
962  : patchWriter()
963  );
964 
965  forAll(pbm, patchi)
966  {
967  if (isA<cyclicAMIPolyPatch>(pbm[patchi]))
968  {
969  const cyclicAMIPolyPatch& cpp =
970  refCast<const cyclicAMIPolyPatch>(pbm[patchi]);
971 
972  if (cpp.owner())
973  {
974  Info<< "Calculating AMI weights between owner patch: "
975  << cpp.name() << " and neighbour patch: "
976  << cpp.neighbPatch().name() << endl;
977 
978  const AMIPatchToPatchInterpolation& ami =
979  cpp.AMI();
980 
981  {
982  // Collect geometry
983  labelList pointToGlobal;
984  labelList uniqueMeshPointLabels;
985  autoPtr<globalIndex> globalPoints;
986  autoPtr<globalIndex> globalFaces;
987  faceList mergedFaces;
988  pointField mergedPoints;
990  (
991  mesh,
992  cpp.localFaces(),
993  cpp.meshPoints(),
994  cpp.meshPointMap(),
995 
996  pointToGlobal,
997  uniqueMeshPointLabels,
998  globalPoints,
999  globalFaces,
1000 
1001  mergedFaces,
1002  mergedPoints
1003  );
1004  // Collect field
1005  scalarField mergedWeights;
1006  globalFaces().gather
1007  (
1010  ami.srcWeightsSum(),
1011  mergedWeights
1012  );
1013 
1014  if (Pstream::master())
1015  {
1016  wr.write
1017  (
1018  "postProcessing",
1019  "src_" + tmName,
1020  mergedPoints,
1021  mergedFaces,
1022  "weightsSum",
1023  mergedWeights,
1024  false
1025  );
1026  }
1027  }
1028  {
1029  // Collect geometry
1030  labelList pointToGlobal;
1031  labelList uniqueMeshPointLabels;
1032  autoPtr<globalIndex> globalPoints;
1033  autoPtr<globalIndex> globalFaces;
1034  faceList mergedFaces;
1035  pointField mergedPoints;
1037  (
1038  mesh,
1039  cpp.neighbPatch().localFaces(),
1040  cpp.neighbPatch().meshPoints(),
1041  cpp.neighbPatch().meshPointMap(),
1042 
1043  pointToGlobal,
1044  uniqueMeshPointLabels,
1045  globalPoints,
1046  globalFaces,
1047 
1048  mergedFaces,
1049  mergedPoints
1050  );
1051  // Collect field
1052  scalarField mergedWeights;
1053  globalFaces().gather
1054  (
1057  ami.tgtWeightsSum(),
1058  mergedWeights
1059  );
1060 
1061  if (Pstream::master())
1062  {
1063  wr.write
1064  (
1065  "postProcessing",
1066  "tgt_" + tmName,
1067  mergedPoints,
1068  mergedFaces,
1069  "weightsSum",
1070  mergedWeights,
1071  false
1072  );
1073  }
1074  }
1075  }
1076  }
1077  }
1078  }
1079 
1080  return noFailedChecks;
1081 }
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< Face, FaceList, PointField, PointType > &p, Field< PointType > &mergedPoints, List< Face > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
static void syncBoundaryFaceList(const polyMesh &, UList< T > &, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronize values on boundary faces only.
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Unit conversion functions.
uint8_t direction
Definition: direction.H:45
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:418
List< face > faceList
Definition: faceListFwd.H:43
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:412
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:275
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
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 dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
label checkGeometry(const polyMesh &mesh, const bool allGeometry, const autoPtr< surfaceWriter > &, const autoPtr< writer< scalar >> &)
const cellShapeList & cells
const pointField & points
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
bool checkWedges(const polyMesh &, const bool report, const Vector< label > &, labelHashSet *)
Check wedge orientation.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
void mergeAndWrite(const polyMesh &mesh, const surfaceWriter &writer, const word &name, const indirectPrimitivePatch setPatch, const fileName &outputDir)
Generate merged surface on master and write. Needs input patch.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
vector point
Point is a vector.
Definition: point.H:41
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)
virtual Ostream & write(const token &)=0
Write next token to stream.
AMIInterpolation< PrimitivePatch< face, SubList, const pointField & >, PrimitivePatch< face, SubList, const pointField & > > AMIPatchToPatchInterpolation
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.
static List< int > & procID(label communicator)
Process ID of given process index.
Definition: UPstream.H:429
Namespace for OpenFOAM.