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 | Copyright (C) 2011-2016 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 {
175  argList::validArgs.append("surfaceFile");
177  (
178  "checkSelfIntersection",
179  "also check for self-intersection"
180  );
182  (
183  "splitNonManifold",
184  "split surface along non-manifold edges"
185  " (default split is fully disconnected)"
186  );
188  (
189  "verbose",
190  "verbose operation"
191  );
193  (
194  "blockMesh",
195  "write vertices/blocks for blockMeshDict"
196  );
197 
198  argList args(argc, argv);
199 
200  const fileName surfFileName = args[1];
201  const bool checkSelfIntersect = args.optionFound("checkSelfIntersection");
202  const bool verbose = args.optionFound("verbose");
203  const bool splitNonManifold = args.optionFound("splitNonManifold");
204 
205  Info<< "Reading surface from " << surfFileName << " ..." << nl << endl;
206 
207 
208  // Read
209  // ~~~~
210 
211  triSurface surf(surfFileName);
212 
213 
214  Info<< "Statistics:" << endl;
215  surf.writeStats(Info);
216  Info<< endl;
217 
218  // write bounding box corners
219  if (args.optionFound("blockMesh"))
220  {
221  pointField cornerPts(boundBox(surf.points(), false).points());
222 
223  Info<< "// blockMeshDict info" << nl << nl;
224 
225  Info<< "vertices\n(" << nl;
226  forAll(cornerPts, ptI)
227  {
228  Info<< " " << cornerPts[ptI] << nl;
229  }
230 
231  // number of divisions needs adjustment later
232  Info<< ");\n" << nl
233  << "blocks\n"
234  << "(\n"
235  << " hex (0 1 2 3 4 5 6 7) (10 10 10) simpleGrading (1 1 1)\n"
236  << ");\n" << nl;
237 
238  Info<< "edges\n();" << nl
239  << "patches\n();" << endl;
240 
241  Info<< nl << "// end blockMeshDict info" << nl << endl;
242  }
243 
244 
245  // Region sizes
246  // ~~~~~~~~~~~~
247 
248  {
249  labelList regionSize(surf.patches().size(), 0);
250 
251  forAll(surf, facei)
252  {
253  label region = surf[facei].region();
254 
255  if (region < 0 || region >= regionSize.size())
256  {
258  << "Triangle " << facei << " vertices " << surf[facei]
259  << " has region " << region << " which is outside the range"
260  << " of regions 0.." << surf.patches().size()-1
261  << endl;
262  }
263  else
264  {
265  regionSize[region]++;
266  }
267  }
268 
269  Info<< "Region\tSize" << nl
270  << "------\t----" << nl;
271  forAll(surf.patches(), patchi)
272  {
273  Info<< surf.patches()[patchi].name() << '\t'
274  << regionSize[patchi] << nl;
275  }
276  Info<< nl << endl;
277  }
278 
279 
280  // Check triangles
281  // ~~~~~~~~~~~~~~~
282 
283  {
284  DynamicList<label> illegalFaces(surf.size()/100 + 1);
285 
286  forAll(surf, facei)
287  {
288  if (!validTri(verbose, surf, facei))
289  {
290  illegalFaces.append(facei);
291  }
292  }
293 
294  if (illegalFaces.size())
295  {
296  Info<< "Surface has " << illegalFaces.size()
297  << " illegal triangles." << endl;
298 
299  OFstream str("illegalFaces");
300  Info<< "Dumping conflicting face labels to " << str.name() << endl
301  << "Paste this into the input for surfaceSubset" << endl;
302  str << illegalFaces;
303  }
304  else
305  {
306  Info<< "Surface has no illegal triangles." << endl;
307  }
308  Info<< endl;
309  }
310 
311 
312 
313  // Triangle quality
314  // ~~~~~~~~~~~~~~~~
315 
316  {
317  scalarField triQ(surf.size(), 0);
318  forAll(surf, facei)
319  {
320  const labelledTri& f = surf[facei];
321 
322  if (f[0] == f[1] || f[0] == f[2] || f[1] == f[2])
323  {
324  //WarningIn(args.executable())
325  // << "Illegal triangle " << facei << " vertices " << f
326  // << " coords " << f.points(surf.points()) << endl;
327  }
328  else
329  {
330  triQ[facei] = triPointRef
331  (
332  surf.points()[f[0]],
333  surf.points()[f[1]],
334  surf.points()[f[2]]
335  ).quality();
336  }
337  }
338 
339  labelList binCount = countBins(0, 1, 20, triQ);
340 
341  Info<< "Triangle quality (equilateral=1, collapsed=0):"
342  << endl;
343 
344 
345  OSstream& os = Info;
346  os.width(4);
347 
348  scalar dist = (1.0 - 0.0)/20.0;
349  scalar min = 0;
350  forAll(binCount, binI)
351  {
352  Info<< " " << min << " .. " << min+dist << " : "
353  << 1.0/surf.size() * binCount[binI]
354  << endl;
355  min += dist;
356  }
357  Info<< endl;
358 
359  label minIndex = findMin(triQ);
360  label maxIndex = findMax(triQ);
361 
362  Info<< " min " << triQ[minIndex] << " for triangle " << minIndex
363  << nl
364  << " max " << triQ[maxIndex] << " for triangle " << maxIndex
365  << nl
366  << endl;
367 
368 
369  if (triQ[minIndex] < SMALL)
370  {
372  << triQ[minIndex] << ". This might give problems in"
373  << " self-intersection testing later on." << endl;
374  }
375 
376  // Dump for subsetting
377  {
378  DynamicList<label> problemFaces(surf.size()/100+1);
379 
380  forAll(triQ, facei)
381  {
382  if (triQ[facei] < 1e-11)
383  {
384  problemFaces.append(facei);
385  }
386  }
387 
388  if (!problemFaces.empty())
389  {
390  OFstream str("badFaces");
391 
392  Info<< "Dumping bad quality faces to " << str.name() << endl
393  << "Paste this into the input for surfaceSubset" << nl
394  << nl << endl;
395 
396  str << problemFaces;
397  }
398  }
399  }
400 
401 
402 
403  // Edges
404  // ~~~~~
405  {
406  const edgeList& edges = surf.edges();
407  const pointField& localPoints = surf.localPoints();
408 
409  scalarField edgeMag(edges.size());
410 
411  forAll(edges, edgeI)
412  {
413  edgeMag[edgeI] = edges[edgeI].mag(localPoints);
414  }
415 
416  label minEdgeI = findMin(edgeMag);
417  label maxEdgeI = findMax(edgeMag);
418 
419  const edge& minE = edges[minEdgeI];
420  const edge& maxE = edges[maxEdgeI];
421 
422 
423  Info<< "Edges:" << nl
424  << " min " << edgeMag[minEdgeI] << " for edge " << minEdgeI
425  << " points " << localPoints[minE[0]] << localPoints[minE[1]]
426  << nl
427  << " max " << edgeMag[maxEdgeI] << " for edge " << maxEdgeI
428  << " points " << localPoints[maxE[0]] << localPoints[maxE[1]]
429  << nl
430  << endl;
431  }
432 
433 
434 
435  // Close points
436  // ~~~~~~~~~~~~
437  {
438  const edgeList& edges = surf.edges();
439  const pointField& localPoints = surf.localPoints();
440 
441  const boundBox bb(localPoints);
442  scalar smallDim = 1e-6 * bb.mag();
443 
444  Info<< "Checking for points less than 1e-6 of bounding box ("
445  << bb.span() << " metre) apart."
446  << endl;
447 
448  // Sort points
449  SortableList<scalar> sortedMag(mag(localPoints));
450 
451  label nClose = 0;
452 
453  for (label i = 1; i < sortedMag.size(); i++)
454  {
455  label ptI = sortedMag.indices()[i];
456 
457  label prevPtI = sortedMag.indices()[i-1];
458 
459  if (mag(localPoints[ptI] - localPoints[prevPtI]) < smallDim)
460  {
461  // Check if neighbours.
462  const labelList& pEdges = surf.pointEdges()[ptI];
463 
464  label edgeI = -1;
465 
466  forAll(pEdges, i)
467  {
468  const edge& e = edges[pEdges[i]];
469 
470  if (e[0] == prevPtI || e[1] == prevPtI)
471  {
472  // point1 and point0 are connected through edge.
473  edgeI = pEdges[i];
474 
475  break;
476  }
477  }
478 
479  nClose++;
480 
481  if (edgeI == -1)
482  {
483  Info<< " close unconnected points "
484  << ptI << ' ' << localPoints[ptI]
485  << " and " << prevPtI << ' '
486  << localPoints[prevPtI]
487  << " distance:"
488  << mag(localPoints[ptI] - localPoints[prevPtI])
489  << endl;
490  }
491  else
492  {
493  Info<< " small edge between points "
494  << ptI << ' ' << localPoints[ptI]
495  << " and " << prevPtI << ' '
496  << localPoints[prevPtI]
497  << " distance:"
498  << mag(localPoints[ptI] - localPoints[prevPtI])
499  << endl;
500  }
501  }
502  }
503 
504  Info<< "Found " << nClose << " nearby points." << nl
505  << endl;
506  }
507 
508 
509 
510  // Check manifold
511  // ~~~~~~~~~~~~~~
512 
513  DynamicList<label> problemFaces(surf.size()/100 + 1);
514 
515  const labelListList& eFaces = surf.edgeFaces();
516 
517  label nSingleEdges = 0;
518  forAll(eFaces, edgeI)
519  {
520  const labelList& myFaces = eFaces[edgeI];
521 
522  if (myFaces.size() == 1)
523  {
524  problemFaces.append(myFaces[0]);
525 
526  nSingleEdges++;
527  }
528  }
529 
530  label nMultEdges = 0;
531  forAll(eFaces, edgeI)
532  {
533  const labelList& myFaces = eFaces[edgeI];
534 
535  if (myFaces.size() > 2)
536  {
537  forAll(myFaces, myFacei)
538  {
539  problemFaces.append(myFaces[myFacei]);
540  }
541 
542  nMultEdges++;
543  }
544  }
545  problemFaces.shrink();
546 
547  if ((nSingleEdges != 0) || (nMultEdges != 0))
548  {
549  Info<< "Surface is not closed since not all edges connected to "
550  << "two faces:" << endl
551  << " connected to one face : " << nSingleEdges << endl
552  << " connected to >2 faces : " << nMultEdges << endl;
553 
554  Info<< "Conflicting face labels:" << problemFaces.size() << endl;
555 
556  OFstream str("problemFaces");
557 
558  Info<< "Dumping conflicting face labels to " << str.name() << endl
559  << "Paste this into the input for surfaceSubset" << endl;
560 
561  str << problemFaces;
562  }
563  else
564  {
565  Info<< "Surface is closed. All edges connected to two faces." << endl;
566  }
567  Info<< endl;
568 
569 
570 
571  // Check singly connected domain
572  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
573 
574  {
575  boolList borderEdge(surf.nEdges(), false);
576  if (splitNonManifold)
577  {
578  const labelListList& eFaces = surf.edgeFaces();
579  forAll(eFaces, edgeI)
580  {
581  if (eFaces[edgeI].size() > 2)
582  {
583  borderEdge[edgeI] = true;
584  }
585  }
586  }
587 
588  labelList faceZone;
589  label numZones = surf.markZones(borderEdge, faceZone);
590 
591  Info<< "Number of unconnected parts : " << numZones << endl;
592 
593  if (numZones > 1)
594  {
595  Info<< "Splitting surface into parts ..." << endl << endl;
596 
597  fileName surfFileNameBase(surfFileName.name());
598  const word fileType = surfFileNameBase.ext();
599  // Strip extension
600  surfFileNameBase = surfFileNameBase.lessExt();
601  // If extension was .gz strip original extension
602  if (fileType == "gz")
603  {
604  surfFileNameBase = surfFileNameBase.lessExt();
605  }
606 
607 
608  {
609  Info<< "Writing zoning to "
610  << fileName
611  (
612  "zone_"
613  + surfFileNameBase
614  + '.'
615  + vtkSurfaceWriter::typeName
616  )
617  << "..." << endl << endl;
618 
619  // Convert data
620  scalarField scalarFaceZone(faceZone.size());
621  forAll(faceZone, i)
622  {
623  scalarFaceZone[i] = faceZone[i];
624  }
625  faceList faces(surf.size());
626  forAll(surf, i)
627  {
628  faces[i] = surf[i].triFaceFace();
629  }
630 
632  (
633  surfFileName.path(),
634  surfFileNameBase,
635  surf.points(),
636  faces,
637  "zone",
638  scalarFaceZone,
639  false // face based data
640  );
641  }
642 
643 
644  for (label zone = 0; zone < numZones; zone++)
645  {
646  boolList includeMap(surf.size(), false);
647 
648  forAll(faceZone, facei)
649  {
650  if (faceZone[facei] == zone)
651  {
652  includeMap[facei] = true;
653  }
654  }
655 
656  labelList pointMap;
658 
659  triSurface subSurf
660  (
661  surf.subsetMesh
662  (
663  includeMap,
664  pointMap,
665  faceMap
666  )
667  );
668 
669  fileName subName(surfFileNameBase + "_" + name(zone) + ".obj");
670 
671  Info<< "writing part " << zone << " size " << subSurf.size()
672  << " to " << subName << endl;
673 
674  subSurf.write(subName);
675  }
676  }
677  }
678 
679 
680 
681  // Check orientation
682  // ~~~~~~~~~~~~~~~~~
683 
684  labelHashSet borderEdge(surf.size()/1000);
685  PatchTools::checkOrientation(surf, false, &borderEdge);
686 
687  //
688  // Colour all faces into zones using borderEdge
689  //
690  labelList normalZone;
691  label numNormalZones = PatchTools::markZones(surf, borderEdge, normalZone);
692 
693  Info<< endl
694  << "Number of zones (connected area with consistent normal) : "
695  << numNormalZones << endl;
696 
697  if (numNormalZones > 1)
698  {
699  Info<< "More than one normal orientation." << endl;
700  }
701  Info<< endl;
702 
703 
704 
705  // Check self-intersection
706  // ~~~~~~~~~~~~~~~~~~~~~~~
707 
708  if (checkSelfIntersect)
709  {
710  Info<< "Checking self-intersection." << endl;
711 
712  triSurfaceSearch querySurf(surf);
713 
714  const indexedOctree<treeDataTriSurface>& tree = querySurf.tree();
715 
716  OBJstream intStream("selfInterPoints.obj");
717 
718  label nInt = 0;
719 
720  forAll(surf.edges(), edgeI)
721  {
722  const edge& e = surf.edges()[edgeI];
723 
724  pointIndexHit hitInfo
725  (
726  tree.findLine
727  (
728  surf.points()[surf.meshPoints()[e[0]]],
729  surf.points()[surf.meshPoints()[e[1]]],
731  (
732  tree,
733  edgeI
734  )
735  )
736  );
737 
738  if (hitInfo.hit())
739  {
740  intStream.write(hitInfo.hitPoint());
741  nInt++;
742  }
743  }
744 
745  if (nInt == 0)
746  {
747  Info<< "Surface is not self-intersecting" << endl;
748  }
749  else
750  {
751  Info<< "Surface is self-intersecting at " << nInt
752  << " locations." << endl;
753  Info<< "Writing intersection points to " << intStream.name()
754  << endl;
755  }
756 
757  //surfaceIntersection inter(querySurf);
758  //
759  //if (inter.cutEdges().empty() && inter.cutPoints().empty())
760  //{
761  // Info<< "Surface is not self-intersecting" << endl;
762  //}
763  //else
764  //{
765  // Info<< "Surface is self-intersecting" << endl;
766  // Info<< "Writing edges of intersection to selfInter.obj" << endl;
767  //
768  // OFstream intStream("selfInter.obj");
769  // forAll(inter.cutPoints(), cutPointi)
770  // {
771  // const point& pt = inter.cutPoints()[cutPointi];
772  //
773  // intStream << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z()
774  // << endl;
775  // }
776  // forAll(inter.cutEdges(), cutEdgeI)
777  // {
778  // const edge& e = inter.cutEdges()[cutEdgeI];
779  //
780  // intStream << "l " << e.start()+1 << ' ' << e.end()+1 << endl;
781  // }
782  //}
783  Info<< endl;
784  }
785 
786 
787  Info<< "\nEnd\n" << endl;
788 
789  return 0;
790 }
791 
792 
793 // ************************************************************************* //
label markZones(const boolList &borderEdge, labelList &faceZone) const
(size and) fills faceZone with zone of face. Zone is area
Definition: triSurface.C:900
Generic output stream.
Definition: OSstream.H:51
#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
A class for handling file names.
Definition: fileName.H:69
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:68
void writeStats(Ostream &) const
Write some statistics.
Definition: triSurface.C:1067
Output to file stream.
Definition: OFstream.H:81
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Field< PointType > & points() const
Return reference to global points.
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
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
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
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.
const labelListList & pointEdges() const
Return point-edge addressing.
Helper class to search on triSurface.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:179
const pointField & points
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
static label markZones(const PrimitivePatch< Face, FaceList, PointField, PointType > &, const BoolListType &borderEdge, labelList &faceZone)
Size and fills faceZone with zone of face.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Base class for zones.
Definition: zone.H:57
A class for handling words, derived from string.
Definition: word.H:59
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:314
static bool checkOrientation(const PrimitivePatch< Face, FaceList, PointField, PointType > &, const bool report=false, labelHashSet *marked=0)
Check for orientation issues.
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:102
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
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).
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
static const char nl
Definition: Ostream.H:262
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: edgeI.H:163
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 labelListList & faceFaces() const
Return face-face addressing.
label nEdges() const
Return number of edges in patch.
pointField points(const pointField &) const
Return the points corresponding to this face.
Definition: triFaceI.H:128
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.
triSurface subsetMesh(const boolList &include, labelList &pointMap, labelList &faceMap) const
Return new surface. Returns pointMap, faceMap from.
Definition: triSurface.C:991
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const labelListList & edgeFaces() const
Return edge-face addressing.
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:83
Triangulated surface description with patch information.
Definition: triSurface.H:65
Foam::argList args(argc, argv)
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:249
virtual int width() const
Get width of output field.
Definition: OSstream.C:273
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Namespace for OpenFOAM.