surfaceCheck.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  surfaceCheck
26 
27 Description
28  Checks geometric and topological quality of a surface.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "triangle.H"
33 #include "triSurface.H"
34 #include "triSurfaceSearch.H"
35 #include "argList.H"
36 #include "OFstream.H"
37 #include "OBJstream.H"
38 #include "SortableList.H"
39 #include "PatchTools.H"
40 #include "vtkSurfaceWriter.H"
41 
42 using namespace Foam;
43 
44 // Does face use valid vertices?
45 bool validTri
46 (
47  const bool verbose,
48  const triSurface& surf,
49  const label facei
50 )
51 {
52  // Simple check on indices ok.
53 
54  const labelledTri& f = surf[facei];
55 
56  forAll(f, fp)
57  {
58  if (f[fp] < 0 || f[fp] >= surf.points().size())
59  {
61  << "triangle " << facei << " vertices " << f
62  << " uses point indices outside point range 0.."
63  << surf.points().size()-1 << endl;
64  return false;
65  }
66  }
67 
68  if ((f[0] == f[1]) || (f[0] == f[2]) || (f[1] == f[2]))
69  {
71  << "triangle " << facei
72  << " uses non-unique vertices " << f
73  << " coords:" << f.points(surf.points())
74  << endl;
75  return false;
76  }
77 
78  // duplicate triangle check
79 
80  const labelList& fFaces = surf.faceFaces()[facei];
81 
82  // Check if faceNeighbours use same points as this face.
83  // Note: discards normal information - sides of baffle are merged.
84  forAll(fFaces, i)
85  {
86  label nbrFacei = fFaces[i];
87 
88  if (nbrFacei <= facei)
89  {
90  // lower numbered faces already checked
91  continue;
92  }
93 
94  const labelledTri& nbrF = surf[nbrFacei];
95 
96  if
97  (
98  ((f[0] == nbrF[0]) || (f[0] == nbrF[1]) || (f[0] == nbrF[2]))
99  && ((f[1] == nbrF[0]) || (f[1] == nbrF[1]) || (f[1] == nbrF[2]))
100  && ((f[2] == nbrF[0]) || (f[2] == nbrF[1]) || (f[2] == nbrF[2]))
101  )
102  {
104  << "triangle " << facei << " vertices " << f
105  << " has the same vertices as triangle " << nbrFacei
106  << " vertices " << nbrF
107  << " coords:" << f.points(surf.points())
108  << endl;
109 
110  return false;
111  }
112  }
113  return true;
114 }
115 
116 
117 labelList countBins
118 (
119  const scalar min,
120  const scalar max,
121  const label nBins,
122  const scalarField& vals
123 )
124 {
125  scalar dist = nBins/(max - min);
126 
127  labelList binCount(nBins, 0);
128 
129  forAll(vals, i)
130  {
131  scalar val = vals[i];
132 
133  label index = -1;
134 
135  if (Foam::mag(val - min) < small)
136  {
137  index = 0;
138  }
139  else if (val >= max - small)
140  {
141  index = nBins - 1;
142  }
143  else
144  {
145  index = label((val - min)*dist);
146 
147  if ((index < 0) || (index >= nBins))
148  {
150  << "value " << val << " at index " << i
151  << " outside range " << min << " .. " << max << endl;
152 
153  if (index < 0)
154  {
155  index = 0;
156  }
157  else
158  {
159  index = nBins - 1;
160  }
161  }
162  }
163  binCount[index]++;
164  }
165 
166  return binCount;
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 int main(int argc, char *argv[])
173 {
174  #include "removeCaseOptions.H"
175 
176  argList::validArgs.append("surface file");
178  (
179  "checkSelfIntersection",
180  "also check for self-intersection"
181  );
183  (
184  "splitNonManifold",
185  "split surface along non-manifold edges"
186  " (default split is fully disconnected)"
187  );
189  (
190  "verbose",
191  "verbose operation"
192  );
194  (
195  "blockMesh",
196  "write vertices/blocks for blockMeshDict"
197  );
198 
199  argList args(argc, argv);
200 
201  const fileName surfFileName = args[1];
202  const bool checkSelfIntersect = args.optionFound("checkSelfIntersection");
203  const bool verbose = args.optionFound("verbose");
204  const bool splitNonManifold = args.optionFound("splitNonManifold");
205 
206  Info<< "Reading surface from " << surfFileName << " ..." << nl << endl;
207 
208 
209  // Read
210  // ~~~~
211 
212  triSurface surf(surfFileName);
213 
214 
215  Info<< "Statistics:" << endl;
216  surf.writeStats(Info);
217  Info<< endl;
218 
219  // write bounding box corners
220  if (args.optionFound("blockMesh"))
221  {
222  pointField cornerPts(boundBox(surf.points(), false).points());
223 
224  Info<< "// blockMeshDict info" << nl << nl;
225 
226  Info<< "vertices\n(" << nl;
227  forAll(cornerPts, ptI)
228  {
229  Info<< " " << cornerPts[ptI] << nl;
230  }
231 
232  // number of divisions needs adjustment later
233  Info<< ");\n" << nl
234  << "blocks\n"
235  << "(\n"
236  << " hex (0 1 2 3 4 5 6 7) (10 10 10) simpleGrading (1 1 1)\n"
237  << ");\n" << nl;
238 
239  Info<< "edges\n();" << nl
240  << "patches\n();" << endl;
241 
242  Info<< nl << "// end blockMeshDict info" << nl << endl;
243  }
244 
245 
246  // Region sizes
247  // ~~~~~~~~~~~~
248 
249  {
250  labelList regionSize(surf.patches().size(), 0);
251 
252  forAll(surf, facei)
253  {
254  label region = surf[facei].region();
255 
256  if (region < 0 || region >= regionSize.size())
257  {
259  << "Triangle " << facei << " vertices " << surf[facei]
260  << " has region " << region << " which is outside the range"
261  << " of regions 0.." << surf.patches().size()-1
262  << endl;
263  }
264  else
265  {
266  regionSize[region]++;
267  }
268  }
269 
270  Info<< "Region\tSize" << nl
271  << "------\t----" << nl;
272  forAll(surf.patches(), patchi)
273  {
274  Info<< surf.patches()[patchi].name() << '\t'
275  << regionSize[patchi] << nl;
276  }
277  Info<< nl << endl;
278  }
279 
280 
281  // Check triangles
282  // ~~~~~~~~~~~~~~~
283 
284  {
285  DynamicList<label> illegalFaces(surf.size()/100 + 1);
286 
287  forAll(surf, facei)
288  {
289  if (!validTri(verbose, surf, facei))
290  {
291  illegalFaces.append(facei);
292  }
293  }
294 
295  if (illegalFaces.size())
296  {
297  Info<< "Surface has " << illegalFaces.size()
298  << " illegal triangles." << endl;
299 
300  OFstream str("illegalFaces");
301  Info<< "Dumping conflicting face labels to " << str.name() << endl
302  << "Paste this into the input for surfaceSubset" << endl;
303  str << illegalFaces;
304  }
305  else
306  {
307  Info<< "Surface has no illegal triangles." << endl;
308  }
309  Info<< endl;
310  }
311 
312 
313 
314  // Triangle quality
315  // ~~~~~~~~~~~~~~~~
316 
317  {
318  scalarField triQ(surf.size(), 0);
319  forAll(surf, facei)
320  {
321  const labelledTri& f = surf[facei];
322 
323  if (f[0] == f[1] || f[0] == f[2] || f[1] == f[2])
324  {
325  // WarningIn(args.executable())
326  // << "Illegal triangle " << facei << " vertices " << f
327  // << " coords " << f.points(surf.points()) << endl;
328  }
329  else
330  {
331  triQ[facei] = triPointRef
332  (
333  surf.points()[f[0]],
334  surf.points()[f[1]],
335  surf.points()[f[2]]
336  ).quality();
337  }
338  }
339 
340  labelList binCount = countBins(0, 1, 20, triQ);
341 
342  Info<< "Triangle quality (equilateral=1, collapsed=0):"
343  << endl;
344 
345 
346  OSstream& os = Info;
347  os.width(4);
348 
349  scalar dist = (1.0 - 0.0)/20.0;
350  scalar min = 0;
351  forAll(binCount, binI)
352  {
353  Info<< " " << min << " .. " << min+dist << " : "
354  << 1.0/surf.size() * binCount[binI]
355  << endl;
356  min += dist;
357  }
358  Info<< endl;
359 
360  label minIndex = findMin(triQ);
361  label maxIndex = findMax(triQ);
362 
363  Info<< " min " << triQ[minIndex] << " for triangle " << minIndex
364  << nl
365  << " max " << triQ[maxIndex] << " for triangle " << maxIndex
366  << nl
367  << endl;
368 
369 
370  if (triQ[minIndex] < small)
371  {
373  << triQ[minIndex] << ". This might give problems in"
374  << " self-intersection testing later on." << endl;
375  }
376 
377  // Dump for subsetting
378  {
379  DynamicList<label> problemFaces(surf.size()/100+1);
380 
381  forAll(triQ, facei)
382  {
383  if (triQ[facei] < 1e-11)
384  {
385  problemFaces.append(facei);
386  }
387  }
388 
389  if (!problemFaces.empty())
390  {
391  OFstream str("badFaces");
392 
393  Info<< "Dumping bad quality faces to " << str.name() << endl
394  << "Paste this into the input for surfaceSubset" << nl
395  << nl << endl;
396 
397  str << problemFaces;
398  }
399  }
400  }
401 
402 
403 
404  // Edges
405  // ~~~~~
406  {
407  const edgeList& edges = surf.edges();
408  const pointField& localPoints = surf.localPoints();
409 
410  scalarField edgeMag(edges.size());
411 
412  forAll(edges, edgeI)
413  {
414  edgeMag[edgeI] = edges[edgeI].mag(localPoints);
415  }
416 
417  label minEdgeI = findMin(edgeMag);
418  label maxEdgeI = findMax(edgeMag);
419 
420  const edge& minE = edges[minEdgeI];
421  const edge& maxE = edges[maxEdgeI];
422 
423 
424  Info<< "Edges:" << nl
425  << " min " << edgeMag[minEdgeI] << " for edge " << minEdgeI
426  << " points " << localPoints[minE[0]] << localPoints[minE[1]]
427  << nl
428  << " max " << edgeMag[maxEdgeI] << " for edge " << maxEdgeI
429  << " points " << localPoints[maxE[0]] << localPoints[maxE[1]]
430  << nl
431  << endl;
432  }
433 
434 
435 
436  // Close points
437  // ~~~~~~~~~~~~
438  {
439  const edgeList& edges = surf.edges();
440  const pointField& localPoints = surf.localPoints();
441 
442  const boundBox bb(localPoints);
443  scalar smallDim = 1e-6 * bb.mag();
444 
445  Info<< "Checking for points less than 1e-6 of bounding box ("
446  << bb.span() << " metre) apart."
447  << endl;
448 
449  // Sort points
450  SortableList<scalar> sortedMag(mag(localPoints));
451 
452  label nClose = 0;
453 
454  for (label i = 1; i < sortedMag.size(); i++)
455  {
456  label ptI = sortedMag.indices()[i];
457 
458  label prevPtI = sortedMag.indices()[i-1];
459 
460  if (mag(localPoints[ptI] - localPoints[prevPtI]) < smallDim)
461  {
462  // Check if neighbours.
463  const labelList& pEdges = surf.pointEdges()[ptI];
464 
465  label edgeI = -1;
466 
467  forAll(pEdges, i)
468  {
469  const edge& e = edges[pEdges[i]];
470 
471  if (e[0] == prevPtI || e[1] == prevPtI)
472  {
473  // point1 and point0 are connected through edge.
474  edgeI = pEdges[i];
475 
476  break;
477  }
478  }
479 
480  nClose++;
481 
482  if (edgeI == -1)
483  {
484  Info<< " close unconnected points "
485  << ptI << ' ' << localPoints[ptI]
486  << " and " << prevPtI << ' '
487  << localPoints[prevPtI]
488  << " distance:"
489  << mag(localPoints[ptI] - localPoints[prevPtI])
490  << endl;
491  }
492  else
493  {
494  Info<< " small edge between points "
495  << ptI << ' ' << localPoints[ptI]
496  << " and " << prevPtI << ' '
497  << localPoints[prevPtI]
498  << " distance:"
499  << mag(localPoints[ptI] - localPoints[prevPtI])
500  << endl;
501  }
502  }
503  }
504 
505  Info<< "Found " << nClose << " nearby points." << nl
506  << endl;
507  }
508 
509 
510 
511  // Check manifold
512  // ~~~~~~~~~~~~~~
513 
514  DynamicList<label> problemFaces(surf.size()/100 + 1);
515 
516  const labelListList& eFaces = surf.edgeFaces();
517 
518  label nSingleEdges = 0;
519  forAll(eFaces, edgeI)
520  {
521  const labelList& myFaces = eFaces[edgeI];
522 
523  if (myFaces.size() == 1)
524  {
525  problemFaces.append(myFaces[0]);
526 
527  nSingleEdges++;
528  }
529  }
530 
531  label nMultEdges = 0;
532  forAll(eFaces, edgeI)
533  {
534  const labelList& myFaces = eFaces[edgeI];
535 
536  if (myFaces.size() > 2)
537  {
538  forAll(myFaces, myFacei)
539  {
540  problemFaces.append(myFaces[myFacei]);
541  }
542 
543  nMultEdges++;
544  }
545  }
546  problemFaces.shrink();
547 
548  if ((nSingleEdges != 0) || (nMultEdges != 0))
549  {
550  Info<< "Surface is not closed since not all edges connected to "
551  << "two faces:" << endl
552  << " connected to one face : " << nSingleEdges << endl
553  << " connected to >2 faces : " << nMultEdges << endl;
554 
555  Info<< "Conflicting face labels:" << problemFaces.size() << endl;
556 
557  OFstream str("problemFaces");
558 
559  Info<< "Dumping conflicting face labels to " << str.name() << endl
560  << "Paste this into the input for surfaceSubset" << endl;
561 
562  str << problemFaces;
563  }
564  else
565  {
566  Info<< "Surface is closed. All edges connected to two faces." << endl;
567  }
568  Info<< endl;
569 
570 
571 
572  // Check singly connected domain
573  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
574 
575  {
576  boolList borderEdge(surf.nEdges(), false);
577  if (splitNonManifold)
578  {
579  const labelListList& eFaces = surf.edgeFaces();
580  forAll(eFaces, edgeI)
581  {
582  if (eFaces[edgeI].size() > 2)
583  {
584  borderEdge[edgeI] = true;
585  }
586  }
587  }
588 
589  labelList faceZone;
590  label numZones = surf.markZones(borderEdge, faceZone);
591 
592  Info<< "Number of unconnected parts : " << numZones << endl;
593 
594  if (numZones > 1)
595  {
596  Info<< "Splitting surface into parts ..." << endl << endl;
597 
598  fileName surfFileNameBase(surfFileName.name());
599  const word fileType = surfFileNameBase.ext();
600  // Strip extension
601  surfFileNameBase = surfFileNameBase.lessExt();
602  // If extension was .gz strip original extension
603  if (fileType == "gz")
604  {
605  surfFileNameBase = surfFileNameBase.lessExt();
606  }
607 
608 
609  {
610  Info<< "Writing zoning to "
611  << fileName
612  (
613  "zone_"
614  + surfFileNameBase
615  + '.'
616  + vtkSurfaceWriter::typeName
617  )
618  << "..." << endl << endl;
619 
620  // Convert data
621  scalarField scalarFaceZone(faceZone.size());
622  forAll(faceZone, i)
623  {
624  scalarFaceZone[i] = faceZone[i];
625  }
626  faceList faces(surf.size());
627  forAll(surf, i)
628  {
629  faces[i] = surf[i].triFaceFace();
630  }
631 
633  (
634  surfFileName.path(),
635  surfFileNameBase,
636  surf.points(),
637  faces,
638  "zone",
639  scalarFaceZone,
640  false // face based data
641  );
642  }
643 
644 
645  for (label zone = 0; zone < numZones; zone++)
646  {
647  boolList includeMap(surf.size(), false);
648 
649  forAll(faceZone, facei)
650  {
651  if (faceZone[facei] == zone)
652  {
653  includeMap[facei] = true;
654  }
655  }
656 
657  labelList pointMap;
659 
660  triSurface subSurf
661  (
662  surf.subsetMesh
663  (
664  includeMap,
665  pointMap,
666  faceMap
667  )
668  );
669 
670  fileName subName(surfFileNameBase + "_" + name(zone) + ".obj");
671 
672  Info<< "writing part " << zone << " size " << subSurf.size()
673  << " to " << subName << endl;
674 
675  subSurf.write(subName);
676  }
677  }
678  }
679 
680 
681 
682  // Check orientation
683  // ~~~~~~~~~~~~~~~~~
684 
685  labelHashSet borderEdge(surf.size()/1000);
686  PatchTools::checkOrientation(surf, false, &borderEdge);
687 
688  //
689  // Colour all faces into zones using borderEdge
690  //
691  labelList normalZone;
692  label numNormalZones = PatchTools::markZones(surf, borderEdge, normalZone);
693 
694  Info<< endl
695  << "Number of zones (connected area with consistent normal) : "
696  << numNormalZones << endl;
697 
698  if (numNormalZones > 1)
699  {
700  Info<< "More than one normal orientation." << endl;
701  }
702  Info<< endl;
703 
704 
705 
706  // Check self-intersection
707  // ~~~~~~~~~~~~~~~~~~~~~~~
708 
709  if (checkSelfIntersect)
710  {
711  Info<< "Checking self-intersection." << endl;
712 
713  triSurfaceSearch querySurf(surf);
714 
715  const indexedOctree<treeDataTriSurface>& tree = querySurf.tree();
716 
717  OBJstream intStream("selfInterPoints.obj");
718 
719  label nInt = 0;
720 
721  forAll(surf.edges(), edgeI)
722  {
723  const edge& e = surf.edges()[edgeI];
724 
725  pointIndexHit hitInfo
726  (
727  tree.findLine
728  (
729  surf.points()[surf.meshPoints()[e[0]]],
730  surf.points()[surf.meshPoints()[e[1]]],
732  (
733  tree,
734  edgeI
735  )
736  )
737  );
738 
739  if (hitInfo.hit())
740  {
741  intStream.write(hitInfo.hitPoint());
742  nInt++;
743  }
744  }
745 
746  if (nInt == 0)
747  {
748  Info<< "Surface is not self-intersecting" << endl;
749  }
750  else
751  {
752  Info<< "Surface is self-intersecting at " << nInt
753  << " locations." << endl;
754  Info<< "Writing intersection points to " << intStream.name()
755  << endl;
756  }
757 
758  // surfaceIntersection inter(querySurf);
759  //
760  // if (inter.cutEdges().empty() && inter.cutPoints().empty())
761  //{
762  // Info<< "Surface is not self-intersecting" << endl;
763  //}
764  // else
765  //{
766  // Info<< "Surface is self-intersecting" << endl;
767  // Info<< "Writing edges of intersection to selfInter.obj" << endl;
768  //
769  // OFstream intStream("selfInter.obj");
770  // forAll(inter.cutPoints(), cutPointi)
771  // {
772  // const point& pt = inter.cutPoints()[cutPointi];
773  //
774  // intStream << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z()
775  // << endl;
776  // }
777  // forAll(inter.cutEdges(), cutEdgeI)
778  // {
779  // const edge& e = inter.cutEdges()[cutEdgeI];
780  //
781  // intStream << "l " << e.start()+1 << ' ' << e.end()+1 << endl;
782  // }
783  //}
784  Info<< endl;
785  }
786 
787 
788  Info<< "\nEnd\n" << endl;
789 
790  return 0;
791 }
792 
793 
794 // ************************************************************************* //
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.
Generic output stream.
Definition: OSstream.H:51
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
A class for handling file names.
Definition: fileName.H:79
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:80
Output to file stream.
Definition: OFstream.H:82
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
label findMin(const ListType &, const label start=0)
Find index of min element (and less than given element).
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
const Field< PointType > & localPoints() const
Return pointField of points in patch.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
static label markZones(const PrimitivePatch< FaceList, PointField > &, const BoolListType &borderEdge, labelList &faceZone)
Size and fills faceZone with zone of face.
Helper class to search on triSurface.
const labelListList & faceFaces() const
Return face-face addressing.
const pointField & points
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
void writeStats(Ostream &) const
Write some statistics.
Definition: triSurface.C:1352
Base class for zones.
Definition: zone.H:57
A class for handling words, derived from string.
Definition: word.H:59
pointField points(const pointField &) const
Return the points corresponding to this face.
Definition: triFaceI.H:128
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:102
word name() const
Return file name (part beyond last /)
Definition: fileName.C:183
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
const Field< PointType > & points() const
Return reference to global points.
const labelListList & edgeFaces() const
Return edge-face addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: edgeI.H:181
Triangle with additional region number.
Definition: labelledTri.H:57
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
label nEdges() const
Return number of edges in patch.
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
static const char nl
Definition: Ostream.H:265
virtual void write(const fileName &outputDir, const fileName &surfaceName, const pointField &points, const faceList &faces, const bool verbose=false) const
Write single surface geometry to file.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
virtual int width() const
Get width of output field.
Definition: OSstream.C:281
label patchi
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
A surfaceWriter for VTK legacy format.
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:253
virtual Ostream & write(const token &)=0
Write next token to stream.
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:117
label markZones(const boolList &borderEdge, labelList &faceZone) const
(size and) fills faceZone with zone of face. Zone is area
Definition: triSurface.C:998
triSurface subsetMesh(const boolList &include, labelList &pointMap, labelList &faceMap) const
Return new surface. Returns pointMap, faceMap from.
Definition: triSurface.C:1089
Triangulated surface description with patch information.
Definition: triSurface.H:66
Foam::argList args(argc, argv)
Namespace for OpenFOAM.