surfaceFeatures.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) 2018-2024 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 Description
25  Identifies features in a surface geometry and writes them to file,
26  based on control parameters specified by the user.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "argList.H"
31 #include "Time.H"
32 #include "triSurfaceMesh.H"
33 #include "featureEdgeMesh.H"
35 #include "surfaceFeatures.H"
36 #include "triSurfaceFields.H"
37 #include "vtkWritePolyData.H"
38 #include "systemDict.H"
39 
40 using namespace Foam;
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  autoPtr<surfaceFeatures> extractFromFile
47  (
48  const fileName& featureEdgeFile,
49  const triSurface& surf,
50  const Switch& geometricTestOnly
51  )
52  {
53  edgeMesh eMesh(featureEdgeFile);
54 
55  // Sometimes duplicate edges are present. Remove them.
56  eMesh.mergeEdges();
57 
58  Info<< nl << "Reading existing feature edges from file "
59  << featureEdgeFile << nl
60  << "Selecting edges purely based on geometric tests: "
61  << geometricTestOnly.asText() << endl;
62 
64  (
65  new surfaceFeatures
66  (
67  surf,
68  eMesh.points(),
69  eMesh.edges(),
70  1e-6,
71  geometricTestOnly
72  )
73  );
74  }
75 
76 
77  autoPtr<surfaceFeatures> extractFromSurface
78  (
79  const triSurface& surf,
80  const Switch& geometricTestOnly,
81  const scalar includedAngle
82  )
83  {
84  Info<< nl
85  << "Constructing feature set from included angle "
86  << includedAngle << nl
87  << "Selecting edges purely based on geometric tests: "
88  << geometricTestOnly.asText() << endl;
89 
91  (
92  new surfaceFeatures
93  (
94  surf,
95  includedAngle,
96  0,
97  0,
98  geometricTestOnly
99  )
100  );
101  }
102 
103 
104  autoPtr<surfaceFeatures> surfaceFeatureSet
105  (
106  const fileName& surfaceFileName,
107  const triSurface& surf,
108  const dictionary& dict,
109  const scalar includedAngle
110  )
111  {
112  const Switch geometricTestOnly = dict.lookupOrDefault<Switch>
113  (
114  "geometricTestOnly",
115  "no"
116  );
117 
118  if (dict.found("files"))
119  {
120  HashTable<fileName, fileName> fileNames(dict.lookup("files"));
121 
122  if (fileNames.found(surfaceFileName))
123  {
124  return extractFromFile
125  (
126  fileNames[surfaceFileName],
127  surf,
128  geometricTestOnly
129  );
130  }
131  else
132  {
133  return extractFromSurface
134  (
135  surf,
136  geometricTestOnly,
137  includedAngle
138  );
139  }
140  }
141  else
142  {
143  return extractFromSurface
144  (
145  surf,
146  geometricTestOnly,
147  includedAngle
148  );
149  }
150  }
151 
152 
153  void extractFeatures
154  (
155  const fileName& surfaceFileName,
156  const Time& runTime,
157  const dictionary& dict
158  )
159  {
160  const fileName sFeatFileName = surfaceFileName.lessExt().name();
161 
162  Info<< "Surface : " << surfaceFileName << nl << endl;
163 
164  const Switch writeVTK =
165  dict.lookupOrDefault<Switch>("writeVTK", "off");
166  const Switch writeObj =
167  dict.lookupOrDefault<Switch>("writeObj", "off");
168  const Switch verboseObj =
169  dict.lookupOrDefault<Switch>("verboseObj", "off");
170 
171  const Switch curvature =
172  dict.lookupOrDefault<Switch>("curvature", "off");
173  const Switch featureProximity =
174  dict.lookupOrDefault<Switch>("featureProximity", "off");
175 
176 
177  Info<< nl << "Feature line extraction is only valid on closed manifold "
178  << "surfaces." << endl;
179 
180  // Read
181  // ~~~~
182 
183  triSurface surf
184  (
185  runTime.path()
186  /runTime.constant()
188  /surfaceFileName
189  );
190 
191  Info<< "Statistics:" << endl;
192  surf.writeStats(Info);
193  Info<< endl;
194 
195  const faceList faces(surf.faces());
196 
197  // Either construct features from surface & featureAngle or read set.
198  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
199 
200  const scalar includedAngle =
201  dict.lookup<scalar>("includedAngle", unitDegrees);
202 
204  (
205  surfaceFeatureSet
206  (
207  surfaceFileName,
208  surf,
209  dict,
210  includedAngle
211  )
212  );
213 
214  // Trim set
215  // ~~~~~~~~
216 
217  if (dict.isDict("trimFeatures"))
218  {
219  dictionary trimDict = dict.subDict("trimFeatures");
220 
221  scalar minLen =
222  trimDict.lookupOrAddDefault<scalar>("minLen", -great);
223 
224  label minElem = trimDict.lookupOrAddDefault<label>("minElem", 0);
225 
226  // Trim away small groups of features
227  if (minElem > 0 || minLen > 0)
228  {
229  Info<< "Removing features of length < "
230  << minLen << endl;
231  Info<< "Removing features with number of edges < "
232  << minElem << endl;
233 
234  set().trimFeatures(minLen, minElem, includedAngle);
235  }
236  }
237 
238 
239  // Subset
240  // ~~~~~~
241 
242  // Convert to marked edges, points
243  List<surfaceFeatures::edgeStatus> edgeStat(set().toStatus());
244 
245  if (dict.isDict("subsetFeatures"))
246  {
247  const dictionary& subsetDict = dict.subDict
248  (
249  "subsetFeatures"
250  );
251 
252  if (subsetDict.found("insideBox"))
253  {
254  treeBoundBox bb(subsetDict.lookup("insideBox")());
255 
256  Info<< "Selecting edges inside bb " << bb;
257  if (writeObj)
258  {
259  Info << " see insideBox.obj";
260  bb.writeOBJ("insideBox.obj");
261  }
262  Info<< endl;
263 
264  selectBox(surf, bb, true, edgeStat);
265  }
266  else if (subsetDict.found("outsideBox"))
267  {
268  treeBoundBox bb(subsetDict.lookup("outsideBox")());
269 
270  Info<< "Removing all edges inside bb " << bb;
271  if (writeObj)
272  {
273  Info<< " see outsideBox.obj" << endl;
274  bb.writeOBJ("outsideBox.obj");
275  }
276  Info<< endl;
277 
278  selectBox(surf, bb, false, edgeStat);
279  }
280 
281  const Switch nonManifoldEdges =
282  subsetDict.lookupOrDefault<Switch>("nonManifoldEdges", "yes");
283 
284  if (!nonManifoldEdges)
285  {
286  Info<< "Removing all non-manifold edges"
287  << " (edges with > 2 connected faces) unless they"
288  << " cross multiple regions" << endl;
289 
290  selectManifoldEdges(surf, 1e-5, includedAngle, edgeStat);
291  }
292 
293  const Switch openEdges =
294  subsetDict.lookupOrDefault<Switch>("openEdges", "yes");
295 
296  if (!openEdges)
297  {
298  Info<< "Removing all open edges"
299  << " (edges with 1 connected face)" << endl;
300 
301  forAll(edgeStat, edgei)
302  {
303  if (surf.edgeFaces()[edgei].size() == 1)
304  {
305  edgeStat[edgei] = surfaceFeatures::NONE;
306  }
307  }
308  }
309 
310  if (subsetDict.found("plane"))
311  {
312  const plane cutPlane(subsetDict.subDict("plane"));
313 
314  selectCutEdges(surf, cutPlane, edgeStat);
315 
316  Info<< "Only edges that intersect the plane with normal "
317  << cutPlane.normal()
318  << " and base point " << cutPlane.refPoint()
319  << " will be included as feature edges."<< endl;
320  }
321  }
322 
323 
324  surfaceFeatures newSet(surf);
325  newSet.setFromStatus(edgeStat, includedAngle);
326 
327  Info<< nl
328  << "Initial feature set:" << nl
329  << " feature points : " << newSet.featurePoints().size() << nl
330  << " feature edges : " << newSet.featureEdges().size() << nl
331  << " of which" << nl
332  << " region edges : " << newSet.nRegionEdges() << nl
333  << " external edges : " << newSet.nExternalEdges() << nl
334  << " internal edges : " << newSet.nInternalEdges() << nl
335  << endl;
336 
337  boolList surfBaffleRegions(surf.patches().size(), false);
338 
339  wordList surfBaffleNames;
340  dict.readIfPresent("baffles", surfBaffleNames);
341 
342  forAll(surf.patches(), pI)
343  {
344  const word& name = surf.patches()[pI].name();
345 
346  if (findIndex(surfBaffleNames, name) != -1)
347  {
348  Info<< "Adding baffle region " << name << endl;
349  surfBaffleRegions[pI] = true;
350  }
351  }
352 
353  // Extracting and writing a extendedFeatureEdgeMesh
355  (
356  newSet,
357  runTime,
358  sFeatFileName + ".extendedFeatureEdgeMesh",
359  surfBaffleRegions
360  );
361 
362 
363  if (dict.isDict("addFeatures"))
364  {
365  const word addFeName = dict.subDict("addFeatures")["name"];
366  Info<< "Adding (without merging) features from " << addFeName
367  << nl << endl;
368 
369  extendedFeatureEdgeMesh addFeMesh
370  (
371  IOobject
372  (
373  addFeName,
374  runTime.time().constant(),
375  "extendedFeatureEdgeMesh",
376  runTime.time(),
379  )
380  );
381  Info<< "Read " << addFeMesh.name() << nl;
382  addFeMesh.writeStats(Info);
383 
384  feMesh.add(addFeMesh);
385  }
386 
387 
388  Info<< nl
389  << "Final feature set:" << nl;
390  feMesh.writeStats(Info);
391 
392  Info<< nl << "Writing extendedFeatureEdgeMesh to "
393  << feMesh.relativeObjectPath() << endl;
394 
395  mkDir(feMesh.path());
396 
397  if (writeObj)
398  {
399  feMesh.writeObj
400  (
401  feMesh.path()/surfaceFileName.lessExt().name(),
402  verboseObj
403  );
404  }
405 
406  feMesh.write();
407 
408  // Write a featureEdgeMesh for backwards compatibility
409  featureEdgeMesh bfeMesh
410  (
411  IOobject
412  (
413  surfaceFileName.lessExt().name() + ".eMesh",
414  runTime.constant(),
416  runTime,
419  false
420  ),
421  feMesh.points(),
422  feMesh.edges()
423  );
424 
425  Info<< nl << "Writing featureEdgeMesh to "
426  << bfeMesh.relativeObjectPath() << endl;
427 
428  bfeMesh.regIOobject::write();
429 
430  // Data to write out in VTK format
431  wordList writeVTKFieldNames;
432  boolList writeVTKFieldIsPointValues;
433  #define DeclareWriteVTKFieldTypeValues(Type, nullArg) \
434  PtrList<const Field<Type>> writeVTKField##Type##Values;
435  FOR_ALL_FIELD_TYPES(DeclareWriteVTKFieldTypeValues);
436  #undef DeclareWriteVTKFieldTypeValues
437 
438  // Find distance between close features
439  if (dict.isDict("closeness"))
440  {
441  Info<< nl << "Extracting internal and external closeness of "
442  << "surface." << endl;
443 
444  const dictionary& closenessDict = dict.subDict("closeness");
445 
446  const Switch faceCloseness =
447  closenessDict.lookupOrDefault<Switch>("faceCloseness", "off");
448  const Switch pointCloseness =
449  closenessDict.lookupOrDefault<Switch>("pointCloseness", "off");
450 
451  const scalar internalAngleTolerance
452  (
453  closenessDict.lookupOrDefault<scalar>
454  (
455  "internalAngleTolerance",
456  unitDegrees,
457  80
458  )
459  );
460 
461  const scalar externalAngleTolerance
462  (
463  closenessDict.lookupOrDefault<scalar>
464  (
465  "externalAngleTolerance",
466  unitDegrees,
467  80
468  )
469  );
470 
471  // Searchable triSurface
472  const triSurfaceMesh searchSurf
473  (
474  IOobject
475  (
476  sFeatFileName + ".closeness",
477  runTime.constant(),
479  runTime
480  ),
481  surf
482  );
483 
484  if (faceCloseness)
485  {
486  Pair<tmp<triSurfaceScalarField>> closenessFields
487  (
488  searchSurf.extractCloseness
489  (
490  internalAngleTolerance,
491  externalAngleTolerance
492  )
493  );
494 
495  Info<< " writing "
496  << closenessFields.first()->name() << endl;
497  closenessFields.first()->write();
498 
499  Info<< " writing "
500  << closenessFields.second()->name() << endl;
501  closenessFields.second()->write();
502 
503  if (writeVTK)
504  {
505  writeVTKFieldNames.append("internalCloseness");
506  writeVTKFieldIsPointValues.append(false);
507  writeVTKFieldscalarValues.append
508  (
509  new scalarField(closenessFields.first())
510  );
511 
512  writeVTKFieldNames.append("externalCloseness");
513  writeVTKFieldIsPointValues.append(false);
514  writeVTKFieldscalarValues.append
515  (
516  new scalarField(closenessFields.second())
517  );
518  }
519  }
520 
521  if (pointCloseness)
522  {
524  (
525  searchSurf.extractPointCloseness
526  (
527  internalAngleTolerance,
528  externalAngleTolerance
529  )
530  );
531 
532  Info<< " writing "
533  << closenessFields.first()->name() << endl;
534  closenessFields.first()->write();
535 
536  Info<< " writing "
537  << closenessFields.second()->name() << endl;
538  closenessFields.second()->write();
539 
540  if (writeVTK)
541  {
542  const faceList faces(searchSurf.faces());
543  const Map<label>& meshPointMap = searchSurf.meshPointMap();
544 
546  internalClosenessPointField = closenessFields.first();
547 
549  externalClosenessPointField = closenessFields.second();
550 
551  scalarField internalCloseness(searchSurf.nPoints(), great);
552  scalarField externalCloseness(searchSurf.nPoints(), great);
553 
554  forAll(meshPointMap, pi)
555  {
556  internalCloseness[pi] =
557  internalClosenessPointField[meshPointMap[pi]];
558 
559  externalCloseness[pi] =
560  externalClosenessPointField[meshPointMap[pi]];
561  }
562 
563  writeVTKFieldNames.append("internalPointCloseness");
564  writeVTKFieldIsPointValues.append(true);
565  writeVTKFieldscalarValues.append
566  (
567  new scalarField(internalCloseness)
568  );
569 
570  writeVTKFieldNames.append("externalPointCloseness");
571  writeVTKFieldIsPointValues.append(true);
572  writeVTKFieldscalarValues.append
573  (
574  new scalarField(externalCloseness)
575  );
576  }
577  }
578  }
579 
580 
581  if (curvature)
582  {
583  Info<< nl << "Extracting curvature of surface at the points."
584  << endl;
585 
587  (
588  IOobject
589  (
590  sFeatFileName + ".curvature",
591  runTime.constant(),
593  runTime
594  ),
595  surf,
596  dimLength,
597  surf.curvature()
598  );
599 
600  k.write();
601 
602  if (writeVTK)
603  {
604  writeVTKFieldNames.append("curvature");
605  writeVTKFieldIsPointValues.append(true);
606  writeVTKFieldscalarValues.append(new scalarField(k));
607  }
608  }
609 
610 
611  if (featureProximity)
612  {
613  Info<< nl << "Extracting proximity of close feature points and "
614  << "edges to the surface" << endl;
615 
616  const scalar searchDistance =
617  dict.lookup<scalar>("maxFeatureProximity");
618 
619  scalarField featureProximity(surf.size(), searchDistance);
620 
621  forAll(surf, fi)
622  {
623  const triPointRef& tri = surf[fi].tri(surf.points());
624 
625  const Tuple2<point, scalar> circle = tri.circumCircle();
626  const point& c = circle.first();
627  const scalar rSqr =
628  min(sqr(4*circle.second()), sqr(searchDistance));
629 
630  pointIndexHitList hitList;
631 
632  feMesh.allNearestFeatureEdges(c, rSqr, hitList);
633  featureProximity[fi] = min
634  (
635  feMesh.minDisconnectedDist(hitList),
636  featureProximity[fi]
637  );
638 
639  feMesh.allNearestFeaturePoints(c, rSqr, hitList);
640  featureProximity[fi] = min
641  (
642  minDist(hitList),
643  featureProximity[fi]
644  );
645  }
646 
647  triSurfaceScalarField featureProximityField
648  (
649  IOobject
650  (
651  sFeatFileName + ".featureProximity",
652  runTime.constant(),
654  runTime,
657  ),
658  surf,
659  dimLength,
660  featureProximity
661  );
662 
663  featureProximityField.write();
664 
665  if (writeVTK)
666  {
667  writeVTKFieldNames.append("featureProximity");
668  writeVTKFieldIsPointValues.append(false);
669  writeVTKFieldscalarValues.append
670  (
671  new scalarField(featureProximity)
672  );
673  }
674  }
675 
676  if (writeVTK)
677  {
678  #define WriteVTKResizeFieldTypeValues(Type, nullArg) \
679  writeVTKField##Type##Values.resize(writeVTKFieldNames.size());
680  FOR_ALL_FIELD_TYPES(WriteVTKResizeFieldTypeValues)
681  #undef WriteVTKResizeFieldTypeValues
682 
684  (
685  runTime.path()
686  /runTime.constant()
688  /sFeatFileName + "Features.vtk",
689  sFeatFileName,
690  runTime.writeFormat() == IOstream::BINARY,
691  surf.points(),
692  labelList(),
693  labelListList(),
694  faces,
695  writeVTKFieldNames,
696  writeVTKFieldIsPointValues,
697  UPtrList<const Field<label>>(writeVTKFieldNames.size())
698  #define WriteVTKFieldTypeValuesParameter(Type, nullArg) \
699  , UPtrList<const Field<Type>>(writeVTKField##Type##Values)
700  FOR_ALL_FIELD_TYPES(WriteVTKFieldTypeValuesParameter)
701  #undef WriteVTKFieldTypeValuesParameter
702  );
703  }
704 
705  Info<< endl;
706  }
707 
708 
709  void extractFeatures
710  (
711  const fileNameList& surfaceFileNames,
712  const Time& runTime,
713  const dictionary& dict
714  )
715  {
716  forAll(surfaceFileNames, i)
717  {
718  extractFeatures(surfaceFileNames[i], runTime, dict);
719  }
720  }
721 }
722 
723 
724 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
725 
726 int main(int argc, char *argv[])
727 {
729  (
730  "extract and write surface features to file"
731  );
733 
734  #include "addDictOption.H"
735 
736  #include "setRootCase.H"
737  #include "createTime.H"
738 
739  const dictionary dict(systemDict("surfaceFeaturesDict", args, runTime));
740 
741  if (dict.found("surfaces"))
742  {
743  extractFeatures
744  (
745  fileNameList(dict.lookup("surfaces")),
746  runTime,
747  dict
748  );
749  }
750  else
751  {
753  {
754  if (!iter().isDict())
755  {
756  continue;
757  }
758 
759  extractFeatures
760  (
761  fileNameList(iter().dict().lookup("surfaces")),
762  runTime,
763  iter().dict()
764  );
765  }
766  }
767 
768  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
769  << " ClockTime = " << runTime.elapsedClockTime() << " s"
770  << nl << endl;
771 
772  Info<< "End\n" << endl;
773 
774  return 0;
775 }
776 
777 
778 
779 // ************************************************************************* //
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An STL-conforming hash table.
Definition: HashTable.H:127
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual Ostream & write(const char)=0
Write character.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:65
const Field< PointType > & points() const
Return reference to global points.
const labelListList & edgeFaces() const
Return edge-face addressing.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
const char * asText() const
Return a text representation of the Switch.
Definition: Switch.C:129
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
fileName path() const
Explicitly inherit path from TimePaths to disambiguate from.
Definition: TimePaths.H:138
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:287
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
const Type2 & second() const
Return second.
Definition: Tuple2.H:131
const Type1 & first() const
Return first.
Definition: Tuple2.H:119
T & first()
Return the first element of the list.
Definition: UListI.H:114
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
time_t elapsedClockTime() const
Returns wall-clock time from clock instantiation.
Definition: clock.C:122
double elapsedCpuTime() const
Return CPU time (in seconds) from the start.
Definition: cpuTime.C:53
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:797
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
T lookupOrAddDefault(const word &, const T &, bool recursive=false, bool patternMatch=true)
Find and return a T, if not found return the given.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
Points connected by edges.
Definition: edgeMesh.H:72
A class for handling file names.
Definition: fileName.H:82
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
fileName lessExt() const
Return file name without extension (part before last .)
Definition: fileName.C:284
const Time & time() const
Return time.
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:61
static const word & geometryDir()
Return the geometry directory name.
Holds feature edges/points of surface.
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
A surface geometry formed of discrete facets, e.g. triangles and/or quadrilaterals,...
Triangulated surface description with patch information.
Definition: triSurface.H:69
tmp< scalarField > curvature() const
Return the curvature of surface at the points.
Definition: triSurface.C:1144
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
void writeStats(Ostream &) const
Write some statistics.
Definition: triSurface.C:1353
faceList faces() const
Return the list of triangles as a faceList.
Definition: triSurface.C:1131
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:80
Tuple2< Point, scalar > circumCircle() const
Return the circum centre and radius.
Definition: triangleI.H:120
A class for handling words, derived from string.
Definition: word.H:62
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const dimensionedScalar c
Speed of light in a vacuum.
void write(const fileName &file, const word &title, const bool binary, const PointField &points, const VertexList &vertices, const LineList &lines, const FaceList &faces, const wordList &fieldNames, const boolList &fieldIsPointValues, const UPtrList< const Field< label >> &fieldLabelValues #define FieldTypeValuesConstArg(Type, nullArg))
Write VTK polygonal data to a file. Takes a PtrList of fields of labels and.
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
scalar minDist(const List< pointIndexHit > &hitList)
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
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void selectBox(const triSurface &surf, const boundBox &bb, const bool removeInside, List< surfaceFeatures::edgeStatus > &edgeStat)
Select edges inside or outside bounding box.
messageStream Info
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
void selectManifoldEdges(const triSurface &surf, const scalar tol, const scalar includedAngle, List< surfaceFeatures::edgeStatus > &edgeStat)
Select manifold edges.
void writeVTK(OFstream &os, const Type &value)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void selectCutEdges(const triSurface &surf, const plane &cutPlane, List< surfaceFeatures::edgeStatus > &edgeStat)
Select edges that are intersected by the given plane.
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
static const char nl
Definition: Ostream.H:266
const unitConversion unitDegrees
dictionary dict
Foam::argList args(argc, argv)
Fields for triSurface.