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