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 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 "vtkSurfaceWriter.H"
38 #include "IOdictionary.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  const Switch closeness =
176  dict.lookupOrDefault<Switch>("closeness", "off");
177 
178 
179  Info<< nl << "Feature line extraction is only valid on closed manifold "
180  << "surfaces." << endl;
181 
182  // Read
183  // ~~~~
184 
185  triSurface surf(runTime.constantPath()/"triSurface"/surfaceFileName);
186 
187  Info<< "Statistics:" << endl;
188  surf.writeStats(Info);
189  Info<< endl;
190 
191  const faceList faces(surf.faces());
192 
193  // Either construct features from surface & featureAngle or read set.
194  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
195 
196  const scalar includedAngle = readScalar(dict.lookup("includedAngle"));
197 
199  (
200  surfaceFeatureSet
201  (
202  surfaceFileName,
203  surf,
204  dict,
205  includedAngle
206  )
207  );
208 
209  // Trim set
210  // ~~~~~~~~
211 
212  if (dict.isDict("trimFeatures"))
213  {
214  dictionary trimDict = dict.subDict("trimFeatures");
215 
216  scalar minLen =
217  trimDict.lookupOrAddDefault<scalar>("minLen", -great);
218 
219  label minElem = trimDict.lookupOrAddDefault<label>("minElem", 0);
220 
221  // Trim away small groups of features
222  if (minElem > 0 || minLen > 0)
223  {
224  Info<< "Removing features of length < "
225  << minLen << endl;
226  Info<< "Removing features with number of edges < "
227  << minElem << endl;
228 
229  set().trimFeatures(minLen, minElem, includedAngle);
230  }
231  }
232 
233 
234  // Subset
235  // ~~~~~~
236 
237  // Convert to marked edges, points
238  List<surfaceFeatures::edgeStatus> edgeStat(set().toStatus());
239 
240  if (dict.isDict("subsetFeatures"))
241  {
242  const dictionary& subsetDict = dict.subDict
243  (
244  "subsetFeatures"
245  );
246 
247  if (subsetDict.found("insideBox"))
248  {
249  treeBoundBox bb(subsetDict.lookup("insideBox")());
250 
251  Info<< "Selecting edges inside bb " << bb;
252  if (writeObj)
253  {
254  Info << " see insideBox.obj";
255  bb.writeOBJ("insideBox.obj");
256  }
257  Info<< endl;
258 
259  selectBox(surf, bb, true, edgeStat);
260  }
261  else if (subsetDict.found("outsideBox"))
262  {
263  treeBoundBox bb(subsetDict.lookup("outsideBox")());
264 
265  Info<< "Removing all edges inside bb " << bb;
266  if (writeObj)
267  {
268  Info<< " see outsideBox.obj" << endl;
269  bb.writeOBJ("outsideBox.obj");
270  }
271  Info<< endl;
272 
273  selectBox(surf, bb, false, edgeStat);
274  }
275 
276  const Switch nonManifoldEdges =
277  subsetDict.lookupOrDefault<Switch>("nonManifoldEdges", "yes");
278 
279  if (!nonManifoldEdges)
280  {
281  Info<< "Removing all non-manifold edges"
282  << " (edges with > 2 connected faces) unless they"
283  << " cross multiple regions" << endl;
284 
285  selectManifoldEdges(surf, 1e-5, includedAngle, edgeStat);
286  }
287 
288  const Switch openEdges =
289  subsetDict.lookupOrDefault<Switch>("openEdges", "yes");
290 
291  if (!openEdges)
292  {
293  Info<< "Removing all open edges"
294  << " (edges with 1 connected face)" << endl;
295 
296  forAll(edgeStat, edgei)
297  {
298  if (surf.edgeFaces()[edgei].size() == 1)
299  {
300  edgeStat[edgei] = surfaceFeatures::NONE;
301  }
302  }
303  }
304 
305  if (subsetDict.found("plane"))
306  {
307  const plane cutPlane(subsetDict.lookup("plane")());
308 
309  selectCutEdges(surf, cutPlane, edgeStat);
310 
311  Info<< "Only edges that intersect the plane with normal "
312  << cutPlane.normal()
313  << " and base point " << cutPlane.refPoint()
314  << " will be included as feature edges."<< endl;
315  }
316  }
317 
318 
319  surfaceFeatures newSet(surf);
320  newSet.setFromStatus(edgeStat, includedAngle);
321 
322  Info<< nl
323  << "Initial feature set:" << nl
324  << " feature points : " << newSet.featurePoints().size() << nl
325  << " feature edges : " << newSet.featureEdges().size() << nl
326  << " of which" << nl
327  << " region edges : " << newSet.nRegionEdges() << nl
328  << " external edges : " << newSet.nExternalEdges() << nl
329  << " internal edges : " << newSet.nInternalEdges() << nl
330  << endl;
331 
332  boolList surfBaffleRegions(surf.patches().size(), false);
333 
334  wordList surfBaffleNames;
335  dict.readIfPresent("baffles", surfBaffleNames);
336 
337  forAll(surf.patches(), pI)
338  {
339  const word& name = surf.patches()[pI].name();
340 
341  if (findIndex(surfBaffleNames, name) != -1)
342  {
343  Info<< "Adding baffle region " << name << endl;
344  surfBaffleRegions[pI] = true;
345  }
346  }
347 
348  // Extracting and writing a extendedFeatureEdgeMesh
350  (
351  newSet,
352  runTime,
353  sFeatFileName + ".extendedFeatureEdgeMesh",
354  surfBaffleRegions
355  );
356 
357 
358  if (dict.isDict("addFeatures"))
359  {
360  const word addFeName = dict.subDict("addFeatures")["name"];
361  Info<< "Adding (without merging) features from " << addFeName
362  << nl << endl;
363 
364  extendedFeatureEdgeMesh addFeMesh
365  (
366  IOobject
367  (
368  addFeName,
369  runTime.time().constant(),
370  "extendedFeatureEdgeMesh",
371  runTime.time(),
374  )
375  );
376  Info<< "Read " << addFeMesh.name() << nl;
377  addFeMesh.writeStats(Info);
378 
379  feMesh.add(addFeMesh);
380  }
381 
382 
383  Info<< nl
384  << "Final feature set:" << nl;
385  feMesh.writeStats(Info);
386 
387  Info<< nl << "Writing extendedFeatureEdgeMesh to "
388  << feMesh.objectPath() << endl;
389 
390  mkDir(feMesh.path());
391 
392  if (writeObj)
393  {
394  feMesh.writeObj
395  (
396  feMesh.path()/surfaceFileName.lessExt().name(),
397  verboseObj
398  );
399  }
400 
401  feMesh.write();
402 
403  // Write a featureEdgeMesh for backwards compatibility
404  featureEdgeMesh bfeMesh
405  (
406  IOobject
407  (
408  surfaceFileName.lessExt().name() + ".eMesh",
409  runTime.constant(),
410  "triSurface",
411  runTime,
414  false
415  ),
416  feMesh.points(),
417  feMesh.edges()
418  );
419 
420  Info<< nl << "Writing featureEdgeMesh to "
421  << bfeMesh.objectPath() << endl;
422 
423  bfeMesh.regIOobject::write();
424 
425 
426  // Find distance between close features
427  if (closeness)
428  {
429  Info<< nl << "Extracting internal and external closeness of "
430  << "surface." << endl;
431 
432  // Searchable triSurface
433  const triSurfaceMesh searchSurf
434  (
435  IOobject
436  (
437  sFeatFileName + ".closeness",
438  runTime.constant(),
439  "triSurface",
440  runTime
441  ),
442  surf
443  );
444 
445  {
446  Pair<tmp<triSurfaceScalarField>> closenessFields
447  (
448  searchSurf.extractCloseness()
449  );
450 
451  closenessFields.first()->write();
452  closenessFields.second()->write();
453 
454  if (writeVTK)
455  {
456  const faceList faces(searchSurf.faces());
457 
459  (
460  runTime.constantPath()/"triSurface",// outputDir
461  searchSurf.objectRegistry::name(), // surfaceName
462  searchSurf.points(),
463  faces,
464  "internalCloseness", // fieldName
465  closenessFields.first(),
466  false, // isNodeValues
467  true // verbose
468  );
469 
471  (
472  runTime.constantPath()/"triSurface",// outputDir
473  searchSurf.objectRegistry::name(), // surfaceName
474  searchSurf.points(),
475  faces,
476  "externalCloseness", // fieldName
477  closenessFields.second(),
478  false, // isNodeValues
479  true // verbose
480  );
481  }
482  }
483 
484  {
486  (
487  searchSurf.extractPointCloseness()
488  );
489 
490  closenessFields.first()->write();
491  closenessFields.second()->write();
492 
493  if (writeVTK)
494  {
495  const faceList faces(searchSurf.faces());
496  const Map<label>& meshPointMap = searchSurf.meshPointMap();
497 
499  internalClosenessPointField = closenessFields.first();
500 
502  externalClosenessPointField = closenessFields.second();
503 
504  scalarField internalCloseness(searchSurf.nPoints(), great);
505  scalarField externalCloseness(searchSurf.nPoints(), great);
506 
507  forAll(meshPointMap, pi)
508  {
509  internalCloseness[pi] =
510  internalClosenessPointField[meshPointMap[pi]];
511 
512  externalCloseness[pi] =
513  externalClosenessPointField[meshPointMap[pi]];
514  }
515 
517  (
518  runTime.constantPath()/"triSurface",// outputDir
519  searchSurf.objectRegistry::name(), // surfaceName
520  searchSurf.points(),
521  faces,
522  "internalPointCloseness", // fieldName
523  internalCloseness,
524  true, // isNodeValues
525  true // verbose
526  );
527 
529  (
530  runTime.constantPath()/"triSurface",// outputDir
531  searchSurf.objectRegistry::name(), // surfaceName
532  searchSurf.points(),
533  faces,
534  "externalPointCloseness", // fieldName
535  externalCloseness,
536  true, // isNodeValues
537  true // verbose
538  );
539  }
540  }
541  }
542 
543 
544  if (curvature)
545  {
546  Info<< nl << "Extracting curvature of surface at the points."
547  << endl;
548 
550  (
551  IOobject
552  (
553  sFeatFileName + ".curvature",
554  runTime.constant(),
555  "triSurface",
556  runTime
557  ),
558  surf,
559  dimLength,
560  surf.curvature()
561  );
562 
563  k.write();
564 
565  if (writeVTK)
566  {
568  (
569  runTime.constantPath()/"triSurface",// outputDir
570  sFeatFileName, // surfaceName
571  surf.points(),
572  faces,
573  "curvature", // fieldName
574  k,
575  true, // isNodeValues
576  true // verbose
577  );
578  }
579  }
580 
581 
582  if (featureProximity)
583  {
584  Info<< nl << "Extracting proximity of close feature points and "
585  << "edges to the surface" << endl;
586 
587  const scalar searchDistance =
588  readScalar(dict.lookup("maxFeatureProximity"));
589 
590  scalarField featureProximity(surf.size(), searchDistance);
591 
592  forAll(surf, fi)
593  {
594  const triPointRef& tri = surf[fi].tri(surf.points());
595  const point& triCentre = tri.circumCentre();
596 
597  const scalar radiusSqr = min
598  (
599  sqr(4*tri.circumRadius()),
600  sqr(searchDistance)
601  );
602 
603  pointIndexHitList hitList;
604 
605  feMesh.allNearestFeatureEdges(triCentre, radiusSqr, hitList);
606  featureProximity[fi] = min
607  (
608  feMesh.minDisconnectedDist(hitList),
609  featureProximity[fi]
610  );
611 
612  feMesh.allNearestFeaturePoints(triCentre, radiusSqr, hitList);
613  featureProximity[fi] = min
614  (
615  minDist(hitList),
616  featureProximity[fi]
617  );
618  }
619 
620  triSurfaceScalarField featureProximityField
621  (
622  IOobject
623  (
624  sFeatFileName + ".featureProximity",
625  runTime.constant(),
626  "triSurface",
627  runTime,
630  ),
631  surf,
632  dimLength,
633  featureProximity
634  );
635 
636  featureProximityField.write();
637 
638  if (writeVTK)
639  {
641  (
642  runTime.constantPath()/"triSurface",// outputDir
643  sFeatFileName, // surfaceName
644  surf.points(),
645  faces,
646  "featureProximity", // fieldName
647  featureProximity,
648  false, // isNodeValues
649  true // verbose
650  );
651  }
652  }
653 
654  Info<< endl;
655  }
656 
657 
658  void extractFeatures
659  (
660  const fileNameList& surfaceFileNames,
661  const Time& runTime,
662  const dictionary& dict
663  )
664  {
665  forAll(surfaceFileNames, i)
666  {
667  extractFeatures(surfaceFileNames[i], runTime, dict);
668  }
669  }
670 }
671 
672 
673 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
674 
675 int main(int argc, char *argv[])
676 {
678  (
679  "extract and write surface features to file"
680  );
682 
683  #include "addDictOption.H"
684 
685  #include "setRootCase.H"
686  #include "createTime.H"
687 
688  const word dictName("surfaceFeaturesDict");
690 
691  Info<< "Reading " << dictName << nl << endl;
692 
693  const IOdictionary dict(dictIO);
694 
695  if (dict.found("surfaces"))
696  {
697  extractFeatures
698  (
699  fileNameList(dict.lookup("surfaces")),
700  runTime,
701  dict
702  );
703  }
704  else
705  {
706  forAllConstIter(dictionary, dict, iter)
707  {
708  if (!iter().isDict())
709  {
710  continue;
711  }
712 
713  extractFeatures
714  (
715  fileNameList(iter().dict().lookup("surfaces")),
716  runTime,
717  iter().dict()
718  );
719  }
720  }
721 
722  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
723  << " ClockTime = " << runTime.elapsedClockTime() << " s"
724  << nl << endl;
725 
726  Info<< "End\n" << endl;
727 
728  return 0;
729 }
730 
731 
732 
733 // ************************************************************************* //
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:57
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 of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
void selectCutEdges(const triSurface &surf, const plane &cutPlane, List< surfaceFeatures::edgeStatus > &edgeStat)
Select edges that are intersected by the given plane.
T lookupOrAddDefault(const word &, const T &, bool recursive=false, bool patternMatch=true)
Find and return a T, if not found return the given.
Fields for triSurface.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
tmp< scalarField > curvature() const
Return the curvature of surface at the points.
Definition: triSurface.C:1143
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
static void noParallel()
Remove the parallel options.
Definition: argList.C:174
void selectBox(const triSurface &surf, const boundBox &bb, const bool removeInside, List< surfaceFeatures::edgeStatus > &edgeStat)
Select edges inside or outside bounding box.
label k
Boltzmann constant.
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
fileName constantPath() const
Return constant path.
Definition: TimePaths.H:146
T & first()
Return the first element of the list.
Definition: UListI.H:114
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
scalar minDist(const List< pointIndexHit > &hitList)
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:653
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:699
IOoject and searching on triSurface.
void selectManifoldEdges(const triSurface &surf, const scalar tol, const scalar includedAngle, List< surfaceFeatures::edgeStatus > &edgeStat)
Select manifold edges.
stressControl lookup("compactNormalStress") >> compactNormalStress
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
void writeStats(Ostream &) const
Write some statistics.
Definition: triSurface.C:1352
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
word name() const
Return file name (part beyond last /)
Definition: fileName.C:183
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
const Field< PointType > & points() const
Return reference to global points.
const labelListList & edgeFaces() const
Return edge-face addressing.
const word dictName("particleTrackDict")
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
scalar circumRadius() const
Return circum-radius.
Definition: triangleI.H:146
An STL-conforming hash table.
Definition: HashTable.H:61
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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.
const Time & time() const
Return time.
Points connected by edges.
Definition: edgeMesh.H:69
const char * asText() const
Return a text representation of the Switch.
Definition: Switch.C:129
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
Point circumCentre() const
Return circum-centre.
Definition: triangleI.H:119
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:67
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
A surfaceWriter for VTK legacy format.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
fileName lessExt() const
Return file name without extension (part before last .)
Definition: fileName.C:272
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
faceList faces() const
Return the list of triangles as a faceList.
Definition: triSurface.C:1130
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:158
Triangulated surface description with patch information.
Definition: triSurface.H:66
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const Type & first() const
Return first.
Definition: Pair.H:87
Holds feature edges/points of surface.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583