surfaceCheck.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  surfaceCheck
29 
30 Group
31  grpSurfaceUtilities
32 
33 Description
34  Check geometric and topological quality of a surface.
35 
36 Usage
37  \b surfaceCheck [OPTION] surfaceFile
38 
39  Options:
40  - \par -checkSelfIntersection
41  Check for self-intersection.
42 
43  - \par -splitNonManifold
44  Split surface along non-manifold edges.
45 
46  - \par -verbose
47  Extra verbosity.
48 
49  - \par -blockMesh
50  Write vertices/blocks for tight-fitting 1 cell blockMeshDict.
51 
52  - \par -outputThreshold <num files>
53  Upper limit on the number of files written.
54  Prevent surfaces with many disconnected parts from writing many files.
55  The default is 10. A value of 0 suppresses file writing.
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #include "triangle.H"
60 #include "edgeHashes.H"
61 #include "triSurface.H"
62 #include "triSurfaceSearch.H"
63 #include "triSurfaceTools.H"
64 #include "argList.H"
65 #include "OFstream.H"
66 #include "OBJstream.H"
67 #include "SortableList.H"
68 #include "PatchTools.H"
69 #include "vtkSurfaceWriter.H"
70 #include "functionObject.H"
71 #include "DynamicField.H"
72 #include "edgeMesh.H"
73 
74 using namespace Foam;
75 
76 labelList countBins
77 (
78  const scalar min,
79  const scalar max,
80  const label nBins,
81  const scalarField& vals
82 )
83 {
84  scalar dist = nBins/(max - min);
85 
86  labelList binCount(nBins, Zero);
87 
88  forAll(vals, i)
89  {
90  scalar val = vals[i];
91 
92  label index = -1;
93 
94  if (Foam::mag(val - min) < SMALL)
95  {
96  index = 0;
97  }
98  else if (val >= max - SMALL)
99  {
100  index = nBins - 1;
101  }
102  else
103  {
104  index = label((val - min)*dist);
105 
106  if ((index < 0) || (index >= nBins))
107  {
109  << "value " << val << " at index " << i
110  << " outside range " << min << " .. " << max << endl;
111 
112  if (index < 0)
113  {
114  index = 0;
115  }
116  else
117  {
118  index = nBins - 1;
119  }
120  }
121  }
122  binCount[index]++;
123  }
124 
125  return binCount;
126 }
127 
128 
129 
130 void writeZoning
131 (
133  const triSurface& surf,
134  const labelList& faceZone,
135  const word& fieldName,
136  const fileName& surfFilePath,
137  const word& surfFileStem
138 )
139 {
140  // Transcribe faces
141  faceList faces;
142  surf.triFaceFaces(faces);
143 
144  writer.open
145  (
146  surf.points(),
147  faces,
148  (surfFilePath / surfFileStem),
149  false // serial - already merged
150  );
151 
152  fileName outputName = writer.write(fieldName, labelField(faceZone));
153 
154  writer.clear();
155 
156  Info<< "Wrote zoning to " << outputName << nl << endl;
157 }
158 
159 
160 void writeParts
161 (
162  const triSurface& surf,
163  const label nFaceZones,
164  const labelList& faceZone,
165  const fileName& surfFilePath,
166  const word& surfFileStem
167 )
168 {
169  for (label zone = 0; zone < nFaceZones; zone++)
170  {
171  boolList includeMap(surf.size(), false);
172 
173  forAll(faceZone, facei)
174  {
175  if (faceZone[facei] == zone)
176  {
177  includeMap[facei] = true;
178  }
179  }
180 
181  triSurface subSurf(surf.subsetMesh(includeMap));
182 
183  fileName subName
184  (
185  surfFilePath
186  / surfFileStem + "_" + name(zone) + ".obj"
187  );
188 
189  Info<< "writing part " << zone << " size " << subSurf.size()
190  << " to " << subName << endl;
191 
192  subSurf.write(subName);
193  }
194 }
195 
196 
197 void syncEdges(const triSurface& p, labelHashSet& markedEdges)
198 {
199  // See comment below about having duplicate edges
200 
201  const edgeList& edges = p.edges();
202  edgeHashSet edgeSet(2*markedEdges.size());
203 
204  for (const label edgei : markedEdges)
205  {
206  edgeSet.insert(edges[edgei]);
207  }
208 
209  forAll(edges, edgei)
210  {
211  if (edgeSet.found(edges[edgei]))
212  {
213  markedEdges.insert(edgei);
214  }
215  }
216 }
217 
218 
219 void syncEdges(const triSurface& p, boolList& isMarkedEdge)
220 {
221  // See comment below about having duplicate edges
222 
223  const edgeList& edges = p.edges();
224 
225  label n = 0;
226  forAll(isMarkedEdge, edgei)
227  {
228  if (isMarkedEdge[edgei])
229  {
230  n++;
231  }
232  }
233 
234  edgeHashSet edgeSet(2*n);
235 
236  forAll(isMarkedEdge, edgei)
237  {
238  if (isMarkedEdge[edgei])
239  {
240  edgeSet.insert(edges[edgei]);
241  }
242  }
243 
244  forAll(edges, edgei)
245  {
246  if (edgeSet.found(edges[edgei]))
247  {
248  isMarkedEdge[edgei] = true;
249  }
250  }
251 }
252 
253 
254 void writeEdgeSet
255 (
256  const word& setName,
257  const triSurface& surf,
258  const labelUList& edgeSet
259 )
260 {
261  // Get compact edge mesh
262  labelList pointToCompact(surf.nPoints(), -1);
263  DynamicField<point> compactPoints(edgeSet.size());
264  DynamicList<edge> compactEdges(edgeSet.size());
265  for (label edgei : edgeSet)
266  {
267  const edge& e = surf.edges()[edgei];
268  edge compactEdge(-1, -1);
269  forAll(e, ep)
270  {
271  label& compacti = pointToCompact[e[ep]];
272  if (compacti == -1)
273  {
274  compacti = compactPoints.size();
275  label pointi = surf.meshPoints()[e[ep]];
276  compactPoints.append(surf.points()[pointi]);
277  }
278  compactEdge[ep] = compacti;
279  }
280  compactEdges.append(compactEdge);
281  }
282 
283  edgeMesh eMesh(std::move(compactPoints), std::move(compactEdges));
284  eMesh.write(setName);
285 }
286 
287 
288 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289 
290 int main(int argc, char *argv[])
291 {
293  (
294  "Check geometric and topological quality of a surface"
295  );
296 
298  argList::addArgument("input", "The input surface file");
299 
301  (
302  "checkSelfIntersection",
303  "Also check for self-intersection"
304  );
306  (
307  "splitNonManifold",
308  "Split surface along non-manifold edges "
309  "(default split is fully disconnected)"
310  );
312  (
313  "Additional verbosity"
314  );
316  (
317  "blockMesh",
318  "Write vertices/blocks for blockMeshDict"
319  );
321  (
322  "outputThreshold",
323  "number",
324  "Upper limit on the number of files written. "
325  "Default is 10, using 0 suppresses file writing."
326  );
328  (
329  "writeSets",
330  "surfaceFormat",
331  "Reconstruct and write problem triangles/edges in selected format"
332  );
333 
334 
335  argList args(argc, argv);
336 
337  const auto surfName = args.get<fileName>(1);
338  const bool checkSelfIntersect = args.found("checkSelfIntersection");
339  const bool splitNonManifold = args.found("splitNonManifold");
340  const label outputThreshold =
341  args.getOrDefault<label>("outputThreshold", 10);
342  const word surfaceFormat = args.getOrDefault<word>("writeSets", "");
343  const bool writeSets = !surfaceFormat.empty();
344 
345  autoPtr<surfaceWriter> surfWriter;
346  word edgeFormat;
347  if (writeSets)
348  {
349  surfWriter = surfaceWriter::New(surfaceFormat);
350 
351  // Option1: hard-coded format
352  edgeFormat = "obj";
355  //edgeFormat = surfaceFormat;
356  //if (!edgeMesh::canWriteType(edgeFormat))
357  //{
358  // edgeFormat = "obj";
359  //}
360  }
361 
362 
363  Info<< "Reading surface from "
364  << args.relativePath(surfName) << " ..." << nl << endl;
365 
366  // Read
367  // ~~~~
368 
369  triSurface surf(surfName);
370 
371 
372  Info<< "Statistics:" << endl;
373  surf.writeStats(Info);
374  Info<< endl;
375 
376 
377  // Split into path and stem (no extension)
378  const fileName surfFilePath(surfName.path());
379  word surfFileStem(surfName.stem());
380 
381  // If originally ".gz", need to strip extension again
382  if (surfName.has_ext("gz"))
383  {
384  surfFileStem.remove_ext();
385  }
386 
387 
388  // write bounding box corners
389  if (args.found("blockMesh"))
390  {
391  pointField cornerPts(boundBox(surf.points(), false).points());
392 
393  Info<< "// blockMeshDict info" << nl << nl;
394 
395  Info<< "vertices\n(" << nl;
396  forAll(cornerPts, pti)
397  {
398  Info<< " " << cornerPts[pti] << nl;
399  }
400 
401  // number of divisions needs adjustment later
402  Info<< ");\n" << nl
403  << "blocks\n"
404  << "(\n"
405  << " hex (0 1 2 3 4 5 6 7) (10 10 10) simpleGrading (1 1 1)\n"
406  << ");\n" << nl;
407 
408  Info<< "edges\n();" << nl
409  << "patches\n();" << endl;
410 
411  Info<< nl << "// end blockMeshDict info" << nl << endl;
412  }
413 
414 
415  // Region sizes
416  // ~~~~~~~~~~~~
417 
418  {
419  labelList regionSize(surf.patches().size(), Zero);
420 
421  forAll(surf, facei)
422  {
423  label region = surf[facei].region();
424 
425  if (region < 0 || region >= regionSize.size())
426  {
428  << "Triangle " << facei << " vertices " << surf[facei]
429  << " has region " << region << " which is outside the range"
430  << " of regions 0.." << surf.patches().size()-1
431  << endl;
432  }
433  else
434  {
435  regionSize[region]++;
436  }
437  }
438 
439  Info<< "Region\tSize" << nl
440  << "------\t----" << nl;
441  forAll(surf.patches(), patchi)
442  {
443  Info<< surf.patches()[patchi].name() << '\t'
444  << regionSize[patchi] << nl;
445  }
446  Info<< nl << endl;
447  }
448 
449 
450  // Check triangles
451  // ~~~~~~~~~~~~~~~
452 
453  {
454  DynamicList<label> illegalFaces(surf.size()/100 + 1);
455 
456  forAll(surf, facei)
457  {
458  // Check silently
459  if (!triSurfaceTools::validTri(surf, facei, false))
460  {
461  illegalFaces.append(facei);
462  }
463  }
464 
465  if (illegalFaces.size())
466  {
467  Info<< "Surface has " << illegalFaces.size()
468  << " illegal triangles." << endl;
469 
470  if (surfWriter)
471  {
472  boolList isIllegalFace(surf.size(), false);
473  UIndirectList<bool>(isIllegalFace, illegalFaces) = true;
474 
475  triSurface subSurf(surf.subsetMesh(isIllegalFace));
476 
477 
478  // Transcribe faces
479  faceList faces;
480  subSurf.triFaceFaces(faces);
481 
482  surfWriter->open
483  (
484  subSurf.points(),
485  faces,
486  (surfFilePath / surfFileStem),
487  false // serial - already merged
488  );
489 
490  fileName outputName = surfWriter->write
491  (
492  "illegal",
493  scalarField(subSurf.size(), Zero)
494  );
495 
496  surfWriter->clear();
497 
498  Info<< "Wrote illegal triangles to "
499  << outputName << nl << endl;
500  }
501  else if (outputThreshold > 0)
502  {
503  forAll(illegalFaces, i)
504  {
505  // Force warning message
506  triSurfaceTools::validTri(surf, illegalFaces[i], true);
507  if (i >= outputThreshold)
508  {
509  Info<< "Suppressing further warning messages."
510  << " Use -outputThreshold to increase number"
511  << " of warnings." << endl;
512  break;
513  }
514  }
515 
516  OFstream str("illegalFaces");
517  Info<< "Dumping conflicting face labels to " << str.name()
518  << endl
519  << "Paste this into the input for surfaceSubset" << endl;
520  str << illegalFaces;
521  }
522  }
523  else
524  {
525  Info<< "Surface has no illegal triangles." << endl;
526  }
527  Info<< endl;
528  }
529 
530 
531 
532  // Triangle quality
533  // ~~~~~~~~~~~~~~~~
534 
535  {
536  scalarField triQ(surf.size(), Zero);
537  forAll(surf, facei)
538  {
539  const labelledTri& f = surf[facei];
540 
541  if (f[0] == f[1] || f[0] == f[2] || f[1] == f[2])
542  {
543  //WarningInFunction
544  // << "Illegal triangle " << facei << " vertices " << f
545  // << " coords " << f.points(surf.points()) << endl;
546  }
547  else
548  {
549  triQ[facei] = triPointRef
550  (
551  surf.points()[f[0]],
552  surf.points()[f[1]],
553  surf.points()[f[2]]
554  ).quality();
555  }
556  }
557 
558  labelList binCount = countBins(0, 1, 20, triQ);
559 
560  Info<< "Triangle quality (equilateral=1, collapsed=0):"
561  << endl;
562 
563 
564  OSstream& os = Info;
565  os.width(4);
566 
567  scalar dist = (1.0 - 0.0)/20.0;
568  scalar min = 0;
569  forAll(binCount, bini)
570  {
571  Info<< " " << min << " .. " << min+dist << " : "
572  << 1.0/surf.size() * binCount[bini]
573  << endl;
574  min += dist;
575  }
576  Info<< endl;
577 
578  labelPair minMaxIds = findMinMax(triQ);
579 
580  const label minIndex = minMaxIds.first();
581  const label maxIndex = minMaxIds.second();
582 
583  Info<< " min " << triQ[minIndex] << " for triangle " << minIndex
584  << nl
585  << " max " << triQ[maxIndex] << " for triangle " << maxIndex
586  << nl
587  << endl;
588 
589 
590  if (triQ[minIndex] < SMALL)
591  {
593  << "Minimum triangle quality is "
594  << triQ[minIndex] << ". This might give problems in"
595  << " self-intersection testing later on." << endl;
596  }
597 
598  // Dump for subsetting
599  if (surfWriter)
600  {
601  // Transcribe faces
602  faceList faces(surf.size());
603  surf.triFaceFaces(faces);
604 
605  surfWriter->open
606  (
607  surf.points(),
608  faces,
609  (surfFilePath / surfFileStem),
610  false // serial - already merged
611  );
612 
613  fileName outputName = surfWriter->write("quality", triQ);
614 
615  surfWriter->clear();
616 
617  Info<< "Wrote triangle-quality to "
618  << outputName << nl << endl;
619  }
620  else if (outputThreshold > 0)
621  {
622  DynamicList<label> problemFaces(surf.size()/100+1);
623 
624  forAll(triQ, facei)
625  {
626  if (triQ[facei] < 1e-11)
627  {
628  problemFaces.append(facei);
629  }
630  }
631 
632  if (!problemFaces.empty())
633  {
634  OFstream str("badFaces");
635 
636  Info<< "Dumping bad quality faces to " << str.name() << endl
637  << "Paste this into the input for surfaceSubset" << nl
638  << nl << endl;
639 
640  str << problemFaces;
641  }
642  }
643  }
644 
645 
646 
647  // Edges
648  // ~~~~~
649  {
650  const edgeList& edges = surf.edges();
651  const pointField& localPoints = surf.localPoints();
652 
653  scalarField edgeMag(edges.size());
654 
655  forAll(edges, edgei)
656  {
657  edgeMag[edgei] = edges[edgei].mag(localPoints);
658  }
659 
660  labelPair minMaxIds = findMinMax(edgeMag);
661 
662  const label minEdgei = minMaxIds.first();
663  const label maxEdgei = minMaxIds.second();
664 
665  const edge& minE = edges[minEdgei];
666  const edge& maxE = edges[maxEdgei];
667 
668 
669  Info<< "Edges:" << nl
670  << " min " << edgeMag[minEdgei] << " for edge " << minEdgei
671  << " points " << localPoints[minE[0]] << localPoints[minE[1]]
672  << nl
673  << " max " << edgeMag[maxEdgei] << " for edge " << maxEdgei
674  << " points " << localPoints[maxE[0]] << localPoints[maxE[1]]
675  << nl
676  << endl;
677  }
678 
679 
680 
681  // Close points
682  // ~~~~~~~~~~~~
683  {
684  const edgeList& edges = surf.edges();
685  const pointField& localPoints = surf.localPoints();
686 
687  const boundBox bb(localPoints);
688  scalar smallDim = 1e-6 * bb.mag();
689 
690  Info<< "Checking for points less than 1e-6 of bounding box ("
691  << bb.span() << " metre) apart."
692  << endl;
693 
694  // Sort points
695  SortableList<scalar> sortedMag(mag(localPoints));
696 
697  label nClose = 0;
698 
699  for (label i = 1; i < sortedMag.size(); i++)
700  {
701  label pti = sortedMag.indices()[i];
702 
703  label prevPti = sortedMag.indices()[i-1];
704 
705  if (mag(localPoints[pti] - localPoints[prevPti]) < smallDim)
706  {
707  // Check if neighbours.
708  const labelList& pEdges = surf.pointEdges()[pti];
709 
710  label edgei = -1;
711 
712  forAll(pEdges, i)
713  {
714  const edge& e = edges[pEdges[i]];
715 
716  if (e[0] == prevPti || e[1] == prevPti)
717  {
718  // point1 and point0 are connected through edge.
719  edgei = pEdges[i];
720 
721  break;
722  }
723  }
724 
725  if (nClose < outputThreshold)
726  {
727  if (edgei == -1)
728  {
729  Info<< " close unconnected points "
730  << pti << ' ' << localPoints[pti]
731  << " and " << prevPti << ' '
732  << localPoints[prevPti]
733  << " distance:"
734  << mag(localPoints[pti] - localPoints[prevPti])
735  << endl;
736  }
737  else
738  {
739  Info<< " small edge between points "
740  << pti << ' ' << localPoints[pti]
741  << " and " << prevPti << ' '
742  << localPoints[prevPti]
743  << " distance:"
744  << mag(localPoints[pti] - localPoints[prevPti])
745  << endl;
746  }
747  }
748  else if (nClose == outputThreshold)
749  {
750  Info<< "Suppressing further warning messages."
751  << " Use -outputThreshold to increase number"
752  << " of warnings." << endl;
753  }
754 
755  nClose++;
756  }
757  }
758 
759  Info<< "Found " << nClose << " nearby points." << nl
760  << endl;
761  }
762 
763 
764 
765  // Check manifold
766  // ~~~~~~~~~~~~~~
767 
768  DynamicList<label> problemFaces(surf.size()/100 + 1);
769  DynamicList<label> openEdges(surf.nEdges()/100 + 1);
770  DynamicList<label> multipleEdges(surf.nEdges()/100 + 1);
771 
772  const labelListList& edgeFaces = surf.edgeFaces();
773 
774  forAll(edgeFaces, edgei)
775  {
776  const labelList& myFaces = edgeFaces[edgei];
777 
778  if (myFaces.size() == 1)
779  {
780  problemFaces.append(myFaces[0]);
781  openEdges.append(edgei);
782  }
783  }
784 
785  forAll(edgeFaces, edgei)
786  {
787  const labelList& myFaces = edgeFaces[edgei];
788 
789  if (myFaces.size() > 2)
790  {
791  forAll(myFaces, myFacei)
792  {
793  problemFaces.append(myFaces[myFacei]);
794  }
795  multipleEdges.append(edgei);
796  }
797  }
798  problemFaces.shrink();
799 
800  if (openEdges.size() || multipleEdges.size())
801  {
802  Info<< "Surface is not closed since not all edges ("
803  << edgeFaces.size() << ") connected to "
804  << "two faces:" << endl
805  << " connected to one face : " << openEdges.size() << endl
806  << " connected to >2 faces : " << multipleEdges.size() << endl;
807 
808  Info<< "Conflicting face labels:" << problemFaces.size() << endl;
809 
810  if (!edgeFormat.empty() && openEdges.size())
811  {
812  const fileName outputName
813  (
814  surfName.lessExt()
815  + "_open."
816  + edgeFormat
817  );
818  Info<< "Writing open edges to "
819  << args.relativePath(outputName) << " ..." << endl;
820  writeEdgeSet(outputName, surf, openEdges);
821  }
822  if (!edgeFormat.empty() && multipleEdges.size())
823  {
824  const fileName outputName
825  (
826  surfName.lessExt()
827  + "_multiply."
828  + edgeFormat
829  );
830  Info<< "Writing multiply connected edges to "
831  << args.relativePath(outputName) << " ..." << endl;
832  writeEdgeSet(outputName, surf, multipleEdges);
833  }
834  }
835  else
836  {
837  Info<< "Surface is closed. All edges connected to two faces." << endl;
838  }
839  Info<< endl;
840 
841 
842 
843  // Check singly connected domain
844  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
845 
846  {
847  boolList borderEdge(surf.nEdges(), false);
848  if (splitNonManifold)
849  {
850  forAll(edgeFaces, edgei)
851  {
852  if (edgeFaces[edgei].size() > 2)
853  {
854  borderEdge[edgei] = true;
855  }
856  }
857  syncEdges(surf, borderEdge);
858  }
859 
861  label numZones = surf.markZones(borderEdge, faceZone);
862 
863  Info<< "Number of unconnected parts : " << numZones << endl;
864 
865  if (numZones > 1 && outputThreshold > 0)
866  {
867  Info<< "Splitting surface into parts ..." << endl << endl;
868 
869  if (!surfWriter)
870  {
871  surfWriter.reset(new surfaceWriters::vtkWriter());
872  }
873 
874  writeZoning
875  (
876  *surfWriter,
877  surf,
878  faceZone,
879  "zone",
880  surfFilePath,
881  surfFileStem
882  );
883 
884  if (numZones > outputThreshold)
885  {
886  Info<< "Limiting number of files to " << outputThreshold
887  << endl;
888  }
889  writeParts
890  (
891  surf,
892  min(outputThreshold, numZones),
893  faceZone,
894  surfFilePath,
895  surfFileStem
896  );
897  }
898  }
899 
900 
901 
902  // Check orientation
903  // ~~~~~~~~~~~~~~~~~
904 
905  labelHashSet borderEdge(surf.size()/1000);
906  PatchTools::checkOrientation(surf, false, &borderEdge);
907 
908  // Bit strange: if a triangle has two same vertices (illegal!) it will
909  // still have three distinct edges (two of which have the same vertices).
910  // In this case the faceEdges addressing is not symmetric, i.e. a
911  // neighbouring, valid, triangle will have correct addressing so 3 distinct
912  // edges so it will miss one of those two identical edges.
913  // - we don't want to fix this in PrimitivePatch since it is too specific
914  // - instead just make sure we mark all identical edges consistently
915  // when we use them for marking.
916 
917  syncEdges(surf, borderEdge);
918 
919  //
920  // Colour all faces into zones using borderEdge
921  //
922  labelList normalZone;
923  label numNormalZones = PatchTools::markZones(surf, borderEdge, normalZone);
924 
925  Info<< endl
926  << "Number of zones (connected area with consistent normal) : "
927  << numNormalZones << endl;
928 
929  if (numNormalZones > 1)
930  {
931  Info<< "More than one normal orientation." << endl;
932 
933  if (outputThreshold > 0)
934  {
935  if (!surfWriter)
936  {
937  surfWriter.reset(new surfaceWriters::vtkWriter());
938  }
939 
940  writeZoning
941  (
942  *surfWriter,
943  surf,
944  normalZone,
945  "normal",
946  surfFilePath,
947  surfFileStem
948  );
949 
950  if (numNormalZones > outputThreshold)
951  {
952  Info<< "Limiting number of files to " << outputThreshold
953  << endl;
954  }
955  writeParts
956  (
957  surf,
958  min(outputThreshold, numNormalZones),
959  normalZone,
960  surfFilePath,
961  surfFileStem + "_normal"
962  );
963  }
964  }
965  Info<< endl;
966 
967 
968 
969  // Check self-intersection
970  // ~~~~~~~~~~~~~~~~~~~~~~~
971 
972  if (checkSelfIntersect)
973  {
974  Info<< "Checking self-intersection." << endl;
975 
976  triSurfaceSearch querySurf(surf);
977 
978  const indexedOctree<treeDataTriSurface>& tree = querySurf.tree();
979 
980  autoPtr<OBJstream> intStreamPtr;
981  if (outputThreshold > 0)
982  {
983  intStreamPtr.reset(new OBJstream("selfInterPoints.obj"));
984  }
985 
986  label nInt = 0;
987 
988  forAll(surf.edges(), edgei)
989  {
990  const edge& e = surf.edges()[edgei];
991  const point& start = surf.points()[surf.meshPoints()[e[0]]];
992  const point& end = surf.points()[surf.meshPoints()[e[1]]];
993 
994  // Exclude hits of connected triangles
996 
997  pointIndexHit hitInfo(tree.findLineAny(start, end, exclOp));
998 
999  if (hitInfo.hit())
1000  {
1001  nInt++;
1002 
1003  if (intStreamPtr)
1004  {
1005  intStreamPtr().write(hitInfo.point());
1006  }
1007 
1008  // Try and find from other side.
1009  pointIndexHit hitInfo2(tree.findLineAny(end, start, exclOp));
1010 
1011  if (hitInfo2.hit() && hitInfo.index() != hitInfo2.index())
1012  {
1013  nInt++;
1014 
1015  if (intStreamPtr)
1016  {
1017  intStreamPtr().write(hitInfo2.point());
1018  }
1019  }
1020  }
1021  }
1022 
1024  //{
1025  // const pointField& localPoints = surf.localPoints();
1026  //
1027  // const boundBox bb(localPoints);
1028  // scalar smallDim = 1e-6 * bb.mag();
1029  // scalar smallDimSqr = Foam::sqr(smallDim);
1030  //
1031  // const pointField& faceCentres = surf.faceCentres();
1032  // forAll(faceCentres, faceI)
1033  // {
1034  // const point& fc = faceCentres[faceI];
1035  // pointIndexHit hitInfo
1036  // (
1037  // tree.findNearest
1038  // (
1039  // fc,
1040  // smallDimSqr,
1041  // findSelfNearOp(tree, faceI)
1042  // )
1043  // );
1044  //
1045  // if (hitInfo.hit() && intStreamPtr)
1046  // {
1047  // intStreamPtr().write(hitInfo.point());
1048  //
1049  // label nearFaceI = hitInfo.index();
1050  // triPointRef nearTri(surf[nearFaceI].tri(surf.points()));
1051  // triStreamPtr().write
1052  // (
1053  // surf[faceI].tri(surf.points()),
1054  // false
1055  // );
1056  // triStreamPtr().write(nearTri, false);
1057  // nInt++;
1058  // }
1059  // }
1060  //}
1061 
1062 
1063  if (nInt == 0)
1064  {
1065  Info<< "Surface is not self-intersecting" << endl;
1066  }
1067  else
1068  {
1069  Info<< "Surface is self-intersecting at " << nInt
1070  << " locations." << endl;
1071 
1072  if (intStreamPtr)
1073  {
1074  Info<< "Writing intersection points to "
1075  << intStreamPtr().name() << endl;
1076  }
1077  }
1078  Info<< endl;
1079  }
1080 
1081 
1082  Info<< "\nEnd\n" << endl;
1083 
1084  return 0;
1085 }
1086 
1087 
1088 // ************************************************************************* //
const labelListList & pointEdges() const
Return point-edge addressing.
static bool checkOrientation(const PrimitivePatch< FaceList, PointField > &, const bool report=false, labelHashSet *marked=0)
Check for orientation issues.
const T & first() const noexcept
Access the first element.
Definition: Pair.H:136
label nPoints() const
Number of points supporting patch faces.
triSurface subsetMesh(const UList< bool > &include, labelList &pointMap, labelList &faceMap) const
Return a new surface subsetted on the selected faces.
Definition: triSurface.C:872
Generic output stream using a standard (STL) stream.
Definition: OSstream.H:50
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:48
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
labelPair findMinMax(const ListType &input, label start=0)
Linear search for the index of the min/max element, similar to std::minmax_element but for lists and ...
A class for handling file names.
Definition: fileName.H:71
A surfaceWriter for VTK legacy (.vtk) or XML (.vtp) format.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const Field< point_type > & localPoints() const
Return pointField of points in patch.
void writeStats(Ostream &os) const
Write some statistics.
Definition: triSurfaceIO.C:348
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:491
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: SortableList.H:56
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual void open(const fileName &outputPath)
Open for output on specified path, using existing surface.
Output to file stream, using an OSstream.
Definition: OFstream.H:49
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:365
static void noParallel()
Remove the parallel options.
Definition: argList.C:551
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:45
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:300
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:227
const labelList & meshPoints() const
Return labelList of mesh points in patch.
static label markZones(const PrimitivePatch< FaceList, PointField > &, const BoolListType &borderEdge, labelList &faceZone)
Size and fills faceZone with zone of face.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
virtual void clear()
Close any open output, remove association with a surface and expire the writer. The parallel flag rem...
Helper class to search on triSurface.
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:510
word outputName("finiteArea-edges.obj")
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:59
Base class for mesh zones.
Definition: zone.H:59
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:118
Tree tree(triangles.begin(), triangles.end())
static void addVerboseOption(const string &usage, bool advanced=false)
Enable a &#39;verbose&#39; bool option, with usage information.
Definition: argList.C:505
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:376
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
const labelListList & edgeFaces() const
Return edge-face addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Dynamically sized Field.
Definition: DynamicField.H:45
A triFace with additional (region) index.
Definition: labelledTri.H:53
const Field< point_type > & points() const noexcept
Return reference to global points.
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:558
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
fileName relativePath(const fileName &input, const bool caseTag=false) const
Return the input relative to the globalPath by stripping off a leading value of the globalPath...
Definition: argListI.H:87
label nEdges() const
Number of edges in patch.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files...
Definition: OBJstream.H:55
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
void triFaceFaces(List< face > &plainFaceList) const
Create a list of faces from the triFaces.
Definition: triSurface.C:716
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:193
OBJstream os(runTime.globalPath()/outputName)
virtual bool open(const fileName &file, bool parallel=UPstream::parRun())
Open file for writing (creates parent directory).
Mesh data needed to do the Finite Area discretisation.
Definition: edgeFaMesh.H:47
labelList f(nPoints)
static bool validTri(const triSurface &, const label facei, const bool verbose=true)
Check single triangle for (topological) validity.
virtual int width() const
Get width of output field.
Definition: OSstream.C:293
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
Non-pointer based hierarchical recursive searching.
#define WarningInFunction
Report a warning using Foam::Warning.
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:57
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:342
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangleFwd.H:39
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
Base class for surface writers.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual fileName write()=0
Write separate surface geometry to file.
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
label markZones(const boolList &borderEdge, labelList &faceZone) const
(size and) fills faceZone with zone of face. Zone is area
Definition: triSurface.C:789
volScalarField & p
const T & second() const noexcept
Access the second element.
Definition: Pair.H:146
Triangulated surface description with patch information.
Definition: triSurface.H:72
Foam::argList args(argc, argv)
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:80
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
bool remove_ext()
Remove extension, return true if string changed.
Definition: stringI.H:93
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157