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-2024 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 
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  );
313  (
314  "blockMesh",
315  "Write vertices/blocks for blockMeshDict"
316  );
318  (
319  "outputThreshold",
320  "number",
321  "Upper limit on the number of files written. "
322  "Default is 10, using 0 suppresses file writing."
323  );
325  (
326  "writeSets",
327  "surfaceFormat",
328  "Reconstruct and write problem triangles/edges in selected format"
329  );
330 
331 
332  argList args(argc, argv);
333 
334  const auto surfName = args.get<fileName>(1);
335  const bool checkSelfIntersect = args.found("checkSelfIntersection");
336  const bool splitNonManifold = args.found("splitNonManifold");
337  const label outputThreshold =
338  args.getOrDefault<label>("outputThreshold", 10);
339  const word surfaceFormat = args.getOrDefault<word>("writeSets", "");
340  const bool writeSets = !surfaceFormat.empty();
341 
342  autoPtr<surfaceWriter> surfWriter;
343  word edgeFormat;
344  if (writeSets)
345  {
346  surfWriter = surfaceWriter::New(surfaceFormat);
347 
348  // Option1: hard-coded format
349  edgeFormat = "obj";
352  //edgeFormat = surfaceFormat;
353  //if (!edgeMesh::canWriteType(edgeFormat))
354  //{
355  // edgeFormat = "obj";
356  //}
357  }
358 
359 
360  Info<< "Reading surface from "
361  << args.relativePath(surfName) << " ..." << nl << endl;
362 
363  // Read
364  // ~~~~
365 
366  triSurface surf(surfName);
367 
368 
369  Info<< "Statistics:" << endl;
370  surf.writeStats(Info);
371  Info<< endl;
372 
373 
374  // Split into path and stem (no extension)
375  const fileName surfFilePath(surfName.path());
376  word surfFileStem(surfName.stem());
377 
378  // If originally ".gz", need to strip extension again
379  if (surfName.has_ext("gz"))
380  {
381  surfFileStem.remove_ext();
382  }
383 
384 
385  // write bounding box corners
386  if (args.found("blockMesh"))
387  {
388  pointField cornerPts(boundBox(surf.points(), false).points());
389 
390  Info<< "// blockMeshDict info" << nl << nl;
391 
392  Info<< "vertices\n(" << nl;
393  forAll(cornerPts, pti)
394  {
395  Info<< " " << cornerPts[pti] << nl;
396  }
397 
398  // number of divisions needs adjustment later
399  Info<< ");\n" << nl
400  << "blocks\n"
401  << "(\n"
402  << " hex (0 1 2 3 4 5 6 7) (10 10 10) simpleGrading (1 1 1)\n"
403  << ");\n" << nl;
404 
405  Info<< "edges\n();" << nl
406  << "patches\n();" << endl;
407 
408  Info<< nl << "// end blockMeshDict info" << nl << endl;
409  }
410 
411 
412  // Region sizes
413  // ~~~~~~~~~~~~
414 
415  {
416  labelList regionSize(surf.patches().size(), Foam::zero{});
417 
418  forAll(surf, facei)
419  {
420  label region = surf[facei].region();
421 
422  if (region < 0 || region >= regionSize.size())
423  {
425  << "Triangle " << facei << " vertices " << surf[facei]
426  << " has region " << region << " which is outside the range"
427  << " of regions 0.." << surf.patches().size()-1
428  << endl;
429  }
430  else
431  {
432  regionSize[region]++;
433  }
434  }
435 
436  Info<< "Index\tRegion\tSize" << nl
437  << "-----\t------\t----" << nl;
438  forAll(surf.patches(), patchi)
439  {
440  Info<< patchi << '\t'
441  << surf.patches()[patchi].name() << '\t'
442  << regionSize[patchi] << nl;
443  }
444  Info<< nl << endl;
445  }
446 
447 
448  // Check triangles
449  // ~~~~~~~~~~~~~~~
450 
451  {
452  DynamicList<label> illegalFaces(surf.size()/100 + 1);
453 
454  forAll(surf, facei)
455  {
456  // Check silently
457  if (!triSurfaceTools::validTri(surf, facei, false))
458  {
459  illegalFaces.append(facei);
460  }
461  }
462 
463  if (illegalFaces.size())
464  {
465  Info<< "Surface has " << illegalFaces.size()
466  << " illegal triangles." << endl;
467 
468  if (surfWriter)
469  {
470  boolList isIllegalFace(surf.size(), false);
471  UIndirectList<bool>(isIllegalFace, illegalFaces) = true;
472 
473  triSurface subSurf(surf.subsetMesh(isIllegalFace));
474 
475 
476  // Transcribe faces
477  faceList faces;
478  subSurf.triFaceFaces(faces);
479 
480  surfWriter->open
481  (
482  subSurf.points(),
483  faces,
484  (surfFilePath / surfFileStem),
485  false // serial - already merged
486  );
487 
488  fileName outputName = surfWriter->write
489  (
490  "illegal",
491  scalarField(subSurf.size(), Zero)
492  );
493 
494  surfWriter->clear();
495 
496  Info<< "Wrote illegal triangles to "
497  << outputName << nl << endl;
498  }
499  else if (outputThreshold > 0)
500  {
501  forAll(illegalFaces, i)
502  {
503  // Force warning message
504  triSurfaceTools::validTri(surf, illegalFaces[i], true);
505  if (i >= outputThreshold)
506  {
507  Info<< "Suppressing further warning messages."
508  << " Use -outputThreshold to increase number"
509  << " of warnings." << endl;
510  break;
511  }
512  }
513 
514  OFstream str("illegalFaces");
515  Info<< "Dumping conflicting face labels to " << str.name()
516  << endl
517  << "Paste this into the input for surfaceSubset" << endl;
518  str << illegalFaces;
519  }
520  }
521  else
522  {
523  Info<< "Surface has no illegal triangles." << endl;
524  }
525  Info<< endl;
526  }
527 
528 
529 
530  // Triangle quality
531  // ~~~~~~~~~~~~~~~~
532 
533  {
534  scalarField triQ(surf.size(), Zero);
535  forAll(surf, facei)
536  {
537  const labelledTri& f = surf[facei];
538 
539  if (f[0] == f[1] || f[0] == f[2] || f[1] == f[2])
540  {
541  //WarningInFunction
542  // << "Illegal triangle " << facei << " vertices " << f
543  // << " coords " << f.points(surf.points()) << endl;
544  }
545  else
546  {
547  triQ[facei] = triPointRef
548  (
549  surf.points()[f[0]],
550  surf.points()[f[1]],
551  surf.points()[f[2]]
552  ).quality();
553  }
554  }
555 
556  labelList binCount = countBins(0, 1, 20, triQ);
557 
558  Info<< "Triangle quality (equilateral=1, collapsed=0):"
559  << endl;
560 
561 
562  OSstream& os = Info;
563  os.width(4);
564 
565  scalar dist = (1.0 - 0.0)/20.0;
566  scalar min = 0;
567  forAll(binCount, bini)
568  {
569  Info<< " " << min << " .. " << min+dist << " : "
570  << 1.0/surf.size() * binCount[bini]
571  << endl;
572  min += dist;
573  }
574  Info<< endl;
575 
576  labelPair minMaxIds = findMinMax(triQ);
577 
578  const label minIndex = minMaxIds.first();
579  const label maxIndex = minMaxIds.second();
580 
581  Info<< " min " << triQ[minIndex] << " for triangle " << minIndex
582  << nl
583  << " max " << triQ[maxIndex] << " for triangle " << maxIndex
584  << nl
585  << endl;
586 
587 
588  if (triQ[minIndex] < SMALL)
589  {
591  << "Minimum triangle quality is "
592  << triQ[minIndex] << ". This might give problems in"
593  << " self-intersection testing later on." << endl;
594  }
595 
596  // Dump for subsetting
597  if (surfWriter)
598  {
599  // Transcribe faces
600  faceList faces(surf.size());
601  surf.triFaceFaces(faces);
602 
603  surfWriter->open
604  (
605  surf.points(),
606  faces,
607  (surfFilePath / surfFileStem),
608  false // serial - already merged
609  );
610 
611  fileName outputName = surfWriter->write("quality", triQ);
612 
613  surfWriter->clear();
614 
615  Info<< "Wrote triangle-quality to "
616  << outputName << nl << endl;
617  }
618  else if (outputThreshold > 0)
619  {
620  DynamicList<label> problemFaces(surf.size()/100+1);
621 
622  forAll(triQ, facei)
623  {
624  if (triQ[facei] < 1e-11)
625  {
626  problemFaces.append(facei);
627  }
628  }
629 
630  if (!problemFaces.empty())
631  {
632  OFstream str("badFaces");
633 
634  Info<< "Dumping bad quality faces to " << str.name() << endl
635  << "Paste this into the input for surfaceSubset" << nl
636  << nl << endl;
637 
638  str << problemFaces;
639  }
640  }
641  }
642 
643 
644 
645  // Edges
646  // ~~~~~
647  {
648  const edgeList& edges = surf.edges();
649  const pointField& localPoints = surf.localPoints();
650 
651  scalarField edgeMag(edges.size());
652 
653  forAll(edges, edgei)
654  {
655  edgeMag[edgei] = edges[edgei].mag(localPoints);
656  }
657 
658  labelPair minMaxIds = findMinMax(edgeMag);
659 
660  const label minEdgei = minMaxIds.first();
661  const label maxEdgei = minMaxIds.second();
662 
663  const edge& minE = edges[minEdgei];
664  const edge& maxE = edges[maxEdgei];
665 
666 
667  Info<< "Edges:" << nl
668  << " min " << edgeMag[minEdgei] << " for edge " << minEdgei
669  << " points " << localPoints[minE[0]] << localPoints[minE[1]]
670  << nl
671  << " max " << edgeMag[maxEdgei] << " for edge " << maxEdgei
672  << " points " << localPoints[maxE[0]] << localPoints[maxE[1]]
673  << nl
674  << endl;
675  }
676 
677 
678 
679  // Close points
680  // ~~~~~~~~~~~~
681  {
682  const edgeList& edges = surf.edges();
683  const pointField& localPoints = surf.localPoints();
684 
685  const boundBox bb(localPoints);
686  scalar smallDim = 1e-6 * bb.mag();
687 
688  Info<< "Checking for points less than 1e-6 of bounding box ("
689  << bb.span() << " metre) apart."
690  << endl;
691 
692  // Sort points
693  SortableList<scalar> sortedMag(mag(localPoints));
694 
695  label nClose = 0;
696 
697  for (label i = 1; i < sortedMag.size(); i++)
698  {
699  label pti = sortedMag.indices()[i];
700 
701  label prevPti = sortedMag.indices()[i-1];
702 
703  if (mag(localPoints[pti] - localPoints[prevPti]) < smallDim)
704  {
705  // Check if neighbours.
706  const labelList& pEdges = surf.pointEdges()[pti];
707 
708  label edgei = -1;
709 
710  forAll(pEdges, i)
711  {
712  const edge& e = edges[pEdges[i]];
713 
714  if (e[0] == prevPti || e[1] == prevPti)
715  {
716  // point1 and point0 are connected through edge.
717  edgei = pEdges[i];
718 
719  break;
720  }
721  }
722 
723  if (nClose < outputThreshold)
724  {
725  if (edgei == -1)
726  {
727  Info<< " close unconnected points "
728  << pti << ' ' << localPoints[pti]
729  << " and " << prevPti << ' '
730  << localPoints[prevPti]
731  << " distance:"
732  << mag(localPoints[pti] - localPoints[prevPti])
733  << endl;
734  }
735  else
736  {
737  Info<< " small edge between points "
738  << pti << ' ' << localPoints[pti]
739  << " and " << prevPti << ' '
740  << localPoints[prevPti]
741  << " distance:"
742  << mag(localPoints[pti] - localPoints[prevPti])
743  << endl;
744  }
745  }
746  else if (nClose == outputThreshold)
747  {
748  Info<< "Suppressing further warning messages."
749  << " Use -outputThreshold to increase number"
750  << " of warnings." << endl;
751  }
752 
753  nClose++;
754  }
755  }
756 
757  Info<< "Found " << nClose << " nearby points." << nl
758  << endl;
759  }
760 
761 
762 
763  // Check manifold
764  // ~~~~~~~~~~~~~~
765 
766  DynamicList<label> problemFaces(surf.size()/100 + 1);
767  DynamicList<label> openEdges(surf.nEdges()/100 + 1);
768  DynamicList<label> multipleEdges(surf.nEdges()/100 + 1);
769 
770  const labelListList& edgeFaces = surf.edgeFaces();
771 
772  forAll(edgeFaces, edgei)
773  {
774  const labelList& myFaces = edgeFaces[edgei];
775 
776  if (myFaces.size() == 1)
777  {
778  problemFaces.append(myFaces[0]);
779  openEdges.append(edgei);
780  }
781  }
782 
783  forAll(edgeFaces, edgei)
784  {
785  const labelList& myFaces = edgeFaces[edgei];
786 
787  if (myFaces.size() > 2)
788  {
789  forAll(myFaces, myFacei)
790  {
791  problemFaces.append(myFaces[myFacei]);
792  }
793  multipleEdges.append(edgei);
794  }
795  }
796  problemFaces.shrink();
797 
798  if (openEdges.size() || multipleEdges.size())
799  {
800  Info<< "Surface is not closed since not all edges ("
801  << edgeFaces.size() << ") connected to "
802  << "two faces:" << endl
803  << " connected to one face : " << openEdges.size() << endl
804  << " connected to >2 faces : " << multipleEdges.size() << endl;
805 
806  Info<< "Conflicting face labels:" << problemFaces.size() << endl;
807 
808  if (!edgeFormat.empty() && openEdges.size())
809  {
810  const fileName outputName
811  (
812  surfName.lessExt()
813  + "_open."
814  + edgeFormat
815  );
816  Info<< "Writing open edges to "
817  << args.relativePath(outputName) << " ..." << endl;
818  writeEdgeSet(outputName, surf, openEdges);
819  }
820  if (!edgeFormat.empty() && multipleEdges.size())
821  {
822  const fileName outputName
823  (
824  surfName.lessExt()
825  + "_multiply."
826  + edgeFormat
827  );
828  Info<< "Writing multiply connected edges to "
829  << args.relativePath(outputName) << " ..." << endl;
830  writeEdgeSet(outputName, surf, multipleEdges);
831  }
832  }
833  else
834  {
835  Info<< "Surface is closed. All edges connected to two faces." << endl;
836  }
837  Info<< endl;
838 
839 
840 
841  // Check singly connected domain
842  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
843 
844  {
845  boolList borderEdge(surf.nEdges(), false);
846  if (splitNonManifold)
847  {
848  forAll(edgeFaces, edgei)
849  {
850  if (edgeFaces[edgei].size() > 2)
851  {
852  borderEdge[edgei] = true;
853  }
854  }
855  syncEdges(surf, borderEdge);
856  }
857 
859  label numZones = surf.markZones(borderEdge, faceZone);
860 
861  Info<< "Number of unconnected parts : " << numZones << endl;
862 
863  if (numZones > 1 && outputThreshold > 0)
864  {
865  Info<< "Splitting surface into parts ..." << endl << endl;
866 
867  if (!surfWriter)
868  {
869  surfWriter.reset(new surfaceWriters::vtkWriter());
870  }
871 
872  writeZoning
873  (
874  *surfWriter,
875  surf,
876  faceZone,
877  "zone",
878  surfFilePath,
879  surfFileStem
880  );
881 
882  if (numZones > outputThreshold)
883  {
884  Info<< "Limiting number of files to " << outputThreshold
885  << endl;
886  }
887  writeParts
888  (
889  surf,
890  min(outputThreshold, numZones),
891  faceZone,
892  surfFilePath,
893  surfFileStem
894  );
895  }
896  }
897 
898 
899 
900  // Check orientation
901  // ~~~~~~~~~~~~~~~~~
902 
903  labelHashSet borderEdge(surf.size()/1000);
904  PatchTools::checkOrientation(surf, false, &borderEdge);
905 
906  // Bit strange: if a triangle has two same vertices (illegal!) it will
907  // still have three distinct edges (two of which have the same vertices).
908  // In this case the faceEdges addressing is not symmetric, i.e. a
909  // neighbouring, valid, triangle will have correct addressing so 3 distinct
910  // edges so it will miss one of those two identical edges.
911  // - we don't want to fix this in PrimitivePatch since it is too specific
912  // - instead just make sure we mark all identical edges consistently
913  // when we use them for marking.
914 
915  syncEdges(surf, borderEdge);
916 
917  //
918  // Colour all faces into zones using borderEdge
919  //
920  labelList normalZone;
921  label numNormalZones = PatchTools::markZones(surf, borderEdge, normalZone);
922 
923  Info<< endl
924  << "Number of zones (connected area with consistent normal) : "
925  << numNormalZones << endl;
926 
927  if (numNormalZones > 1)
928  {
929  Info<< "More than one normal orientation." << endl;
930 
931  if (outputThreshold > 0)
932  {
933  if (!surfWriter)
934  {
935  surfWriter.reset(new surfaceWriters::vtkWriter());
936  }
937 
938  writeZoning
939  (
940  *surfWriter,
941  surf,
942  normalZone,
943  "normal",
944  surfFilePath,
945  surfFileStem
946  );
947 
948  if (numNormalZones > outputThreshold)
949  {
950  Info<< "Limiting number of files to " << outputThreshold
951  << endl;
952  }
953  writeParts
954  (
955  surf,
956  min(outputThreshold, numNormalZones),
957  normalZone,
958  surfFilePath,
959  surfFileStem + "_normal"
960  );
961  }
962  }
963  Info<< endl;
964 
965 
966 
967  // Check self-intersection
968  // ~~~~~~~~~~~~~~~~~~~~~~~
969 
970  if (checkSelfIntersect)
971  {
972  Info<< "Checking self-intersection." << endl;
973 
974  triSurfaceSearch querySurf(surf);
975 
976  const indexedOctree<treeDataTriSurface>& tree = querySurf.tree();
977 
978  autoPtr<OBJstream> intStreamPtr;
979  if (outputThreshold > 0)
980  {
981  intStreamPtr.reset(new OBJstream("selfInterPoints.obj"));
982  }
983 
984  label nInt = 0;
985 
986  forAll(surf.edges(), edgei)
987  {
988  const edge& e = surf.edges()[edgei];
989  const point& start = surf.points()[surf.meshPoints()[e[0]]];
990  const point& end = surf.points()[surf.meshPoints()[e[1]]];
991 
992  // Exclude hits of connected triangles
994 
995  pointIndexHit hitInfo(tree.findLineAny(start, end, exclOp));
996 
997  if (hitInfo.hit())
998  {
999  nInt++;
1000 
1001  if (intStreamPtr)
1002  {
1003  intStreamPtr().write(hitInfo.point());
1004  }
1005 
1006  // Try and find from other side.
1007  pointIndexHit hitInfo2(tree.findLineAny(end, start, exclOp));
1008 
1009  if (hitInfo2.hit() && hitInfo.index() != hitInfo2.index())
1010  {
1011  nInt++;
1012 
1013  if (intStreamPtr)
1014  {
1015  intStreamPtr().write(hitInfo2.point());
1016  }
1017  }
1018  }
1019  }
1020 
1022  //{
1023  // const pointField& localPoints = surf.localPoints();
1024  //
1025  // const boundBox bb(localPoints);
1026  // scalar smallDim = 1e-6 * bb.mag();
1027  // scalar smallDimSqr = Foam::sqr(smallDim);
1028  //
1029  // const pointField& faceCentres = surf.faceCentres();
1030  // forAll(faceCentres, faceI)
1031  // {
1032  // const point& fc = faceCentres[faceI];
1033  // pointIndexHit hitInfo
1034  // (
1035  // tree.findNearest
1036  // (
1037  // fc,
1038  // smallDimSqr,
1039  // findSelfNearOp(tree, faceI)
1040  // )
1041  // );
1042  //
1043  // if (hitInfo.hit() && intStreamPtr)
1044  // {
1045  // intStreamPtr().write(hitInfo.point());
1046  //
1047  // label nearFaceI = hitInfo.index();
1048  // triPointRef nearTri(surf[nearFaceI].tri(surf.points()));
1049  // triStreamPtr().write
1050  // (
1051  // surf[faceI].tri(surf.points()),
1052  // false
1053  // );
1054  // triStreamPtr().write(nearTri, false);
1055  // nInt++;
1056  // }
1057  // }
1058  //}
1059 
1060 
1061  if (nInt == 0)
1062  {
1063  Info<< "Surface is not self-intersecting" << endl;
1064  }
1065  else
1066  {
1067  Info<< "Surface is self-intersecting at " << nInt
1068  << " locations." << endl;
1069 
1070  if (intStreamPtr)
1071  {
1072  Info<< "Writing intersection points to "
1073  << intStreamPtr().name() << endl;
1074  }
1075  }
1076  Info<< endl;
1077  }
1078 
1079 
1080  Info<< "\nEnd\n" << endl;
1081 
1082  return 0;
1083 }
1084 
1085 
1086 // ************************************************************************* //
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:137
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:873
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:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
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:72
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.
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
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:521
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 as an OSstream, normally using std::ofstream for the actual output...
Definition: OFstream.H:71
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:374
static void noParallel()
Remove the parallel options.
Definition: argList.C:584
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
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:232
void write(const word &fieldName, const UList< Type > &field)
Write primitive field of CellData (Poly or Line) or PointData values.
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:421
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...
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:358
Helper class to search on triSurface.
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:509
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 expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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 vertex labels. This can correspond to a directed graph edge or an edge 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 addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:385
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:584
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
static void addVerboseOption(const string &usage="", bool advanced=false)
Enable a &#39;verbose&#39; bool option, with usage information.
Definition: argList.C:520
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
void triFaceFaces(List< face > &plainFaceList) const
Create a list of faces from the triFaces.
Definition: triSurface.C:717
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
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.
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:52
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:351
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangleFwd.H:39
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
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:790
virtual int width() const override
Get width of output field.
Definition: OSstream.C:322
volScalarField & p
const T & second() const noexcept
Access the second element.
Definition: Pair.H:147
Triangulated surface description with patch information.
Definition: triSurface.H:71
Foam::argList args(argc, argv)
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
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:127