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