surfaceFeatureExtract.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  Identifies features in a surface geometry and writes them to file.
26  This utility is deprecated, having been replaced by surfaceFeatures
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  void writeStats(const extendedFeatureEdgeMesh& fem, Ostream& os)
47  {
48  os << " points : " << fem.points().size() << nl
49  << " of which" << nl
50  << " convex : "
51  << fem.concaveStart() << nl
52  << " concave : "
53  << (fem.mixedStart() - fem.concaveStart()) << nl
54  << " mixed : "
55  << (fem.nonFeatureStart() - fem.mixedStart()) << nl
56  << " non-feature : "
57  << (fem.points().size() - fem.nonFeatureStart()) << nl
58  << " edges : " << fem.edges().size() << nl
59  << " of which" << nl
60  << " external edges : "
61  << fem.internalStart() << nl
62  << " internal edges : "
63  << (fem.flatStart() - fem.internalStart()) << nl
64  << " flat edges : "
65  << (fem.openStart() - fem.flatStart()) << nl
66  << " open edges : "
67  << (fem.multipleStart() - fem.openStart()) << nl
68  << " multiply connected : "
69  << (fem.edges().size() - fem.multipleStart()) << endl;
70  }
71 }
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 int main(int argc, char *argv[])
76 {
78  (
79  "extract and write surface features to file"
80  );
82 
83  #include "addDictOption.H"
84 
85  #include "setRootCase.H"
86  #include "createTime.H"
87 
88  const word dictName("surfaceFeatureExtractDict");
90 
91  Info<< "Reading " << dictName << nl << endl;
92 
93  const IOdictionary dict(dictIO);
94 
96  {
97  if (!iter().isDict())
98  {
99  continue;
100  }
101 
102  const dictionary& surfaceDict = iter().dict();
103 
104  if (!surfaceDict.found("extractionMethod"))
105  {
106  continue;
107  }
108 
109  const word extractionMethod = surfaceDict.lookup("extractionMethod");
110 
111  const fileName surfFileName = iter().keyword();
112  const fileName sFeatFileName = surfFileName.lessExt().name();
113 
114  Info<< "Surface : " << surfFileName << nl << endl;
115 
116  const Switch writeVTK =
117  surfaceDict.lookupOrDefault<Switch>("writeVTK", "off");
118  const Switch writeObj =
119  surfaceDict.lookupOrDefault<Switch>("writeObj", "off");
120 
121  const Switch curvature =
122  surfaceDict.lookupOrDefault<Switch>("curvature", "off");
123  const Switch featureProximity =
124  surfaceDict.lookupOrDefault<Switch>("featureProximity", "off");
125  const Switch closeness =
126  surfaceDict.lookupOrDefault<Switch>("closeness", "off");
127 
128 
129  Info<< nl << "Feature line extraction is only valid on closed manifold "
130  << "surfaces." << endl;
131 
132  // Read
133  // ~~~~
134 
135  triSurface surf(runTime.constantPath()/"triSurface"/surfFileName);
136 
137  Info<< "Statistics:" << endl;
138  surf.writeStats(Info);
139  Info<< endl;
140 
141  const faceList faces(surf.faces());
142 
143  // Either construct features from surface & featureAngle or read set.
144  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
145 
146  autoPtr<surfaceFeatures> set;
147 
148  scalar includedAngle = 0.0;
149 
150  if (extractionMethod == "extractFromFile")
151  {
152  const dictionary& extractFromFileDict =
153  surfaceDict.subDict("extractFromFileCoeffs");
154 
155  const fileName featureEdgeFile =
156  extractFromFileDict.lookup("featureEdgeFile");
157 
158  const Switch geometricTestOnly =
159  extractFromFileDict.lookupOrDefault<Switch>
160  (
161  "geometricTestOnly",
162  "no"
163  );
164 
165  edgeMesh eMesh(featureEdgeFile);
166 
167  // Sometimes duplicate edges are present. Remove them.
168  eMesh.mergeEdges();
169 
170  Info<< nl << "Reading existing feature edges from file "
171  << featureEdgeFile << nl
172  << "Selecting edges purely based on geometric tests: "
173  << geometricTestOnly.asText() << endl;
174 
175  set.set
176  (
177  new surfaceFeatures
178  (
179  surf,
180  eMesh.points(),
181  eMesh.edges(),
182  1e-6,
183  geometricTestOnly
184  )
185  );
186  }
187  else if (extractionMethod == "extractFromSurface")
188  {
189  const dictionary& extractFromSurfaceDict =
190  surfaceDict.subDict("extractFromSurfaceCoeffs");
191 
192  includedAngle =
193  readScalar(extractFromSurfaceDict.lookup("includedAngle"));
194 
195  const Switch geometricTestOnly =
196  extractFromSurfaceDict.lookupOrDefault<Switch>
197  (
198  "geometricTestOnly",
199  "no"
200  );
201 
202  Info<< nl
203  << "Constructing feature set from included angle "
204  << includedAngle << nl
205  << "Selecting edges purely based on geometric tests: "
206  << geometricTestOnly.asText() << endl;
207 
208  set.set
209  (
210  new surfaceFeatures
211  (
212  surf,
213  includedAngle,
214  0,
215  0,
216  geometricTestOnly
217  )
218  );
219  }
220  else
221  {
223  << "No initial feature set. Provide either one"
224  << " of extractFromFile (to read existing set)" << nl
225  << " or extractFromSurface (to construct new set from angle)"
226  << exit(FatalError);
227  }
228 
229 
230  // Trim set
231  // ~~~~~~~~
232 
233  if (surfaceDict.isDict("trimFeatures"))
234  {
235  dictionary trimDict = surfaceDict.subDict("trimFeatures");
236 
237  scalar minLen =
238  trimDict.lookupOrAddDefault<scalar>("minLen", -great);
239 
240  label minElem = trimDict.lookupOrAddDefault<label>("minElem", 0);
241 
242  // Trim away small groups of features
243  if (minElem > 0 || minLen > 0)
244  {
245  Info<< "Removing features of length < "
246  << minLen << endl;
247  Info<< "Removing features with number of edges < "
248  << minElem << endl;
249 
250  set().trimFeatures(minLen, minElem, includedAngle);
251  }
252  }
253 
254 
255  // Subset
256  // ~~~~~~
257 
258  // Convert to marked edges, points
259  List<surfaceFeatures::edgeStatus> edgeStat(set().toStatus());
260 
261  if (surfaceDict.isDict("subsetFeatures"))
262  {
263  const dictionary& subsetDict = surfaceDict.subDict
264  (
265  "subsetFeatures"
266  );
267 
268  if (subsetDict.found("insideBox"))
269  {
270  treeBoundBox bb(subsetDict.lookup("insideBox")());
271 
272  Info<< "Selecting edges inside bb " << bb
273  << " see insideBox.obj" << endl;
274  bb.writeOBJ("insideBox.obj");
275  selectBox(surf, bb, true, edgeStat);
276  }
277  else if (subsetDict.found("outsideBox"))
278  {
279  treeBoundBox bb(subsetDict.lookup("outsideBox")());
280 
281  Info<< "Removing all edges inside bb " << bb
282  << " see outsideBox.obj" << endl;
283  bb.writeOBJ("outsideBox.obj");
284  selectBox(surf, bb, false, edgeStat);
285  }
286 
287  const Switch nonManifoldEdges =
288  subsetDict.lookupOrDefault<Switch>("nonManifoldEdges", "yes");
289 
290  if (!nonManifoldEdges)
291  {
292  Info<< "Removing all non-manifold edges"
293  << " (edges with > 2 connected faces) unless they"
294  << " cross multiple regions" << endl;
295 
296  selectManifoldEdges(surf, 1e-5, includedAngle, edgeStat);
297  }
298 
299  const Switch openEdges =
300  subsetDict.lookupOrDefault<Switch>("openEdges", "yes");
301 
302  if (!openEdges)
303  {
304  Info<< "Removing all open edges"
305  << " (edges with 1 connected face)" << endl;
306 
307  forAll(edgeStat, edgei)
308  {
309  if (surf.edgeFaces()[edgei].size() == 1)
310  {
311  edgeStat[edgei] = surfaceFeatures::NONE;
312  }
313  }
314  }
315 
316  if (subsetDict.found("plane"))
317  {
318  const plane cutPlane(subsetDict.lookup("plane")());
319 
320  selectCutEdges(surf, cutPlane, edgeStat);
321 
322  Info<< "Only edges that intersect the plane with normal "
323  << cutPlane.normal()
324  << " and base point " << cutPlane.refPoint()
325  << " will be included as feature edges."<< endl;
326  }
327  }
328 
329 
330  surfaceFeatures newSet(surf);
331  newSet.setFromStatus(edgeStat, includedAngle);
332 
333  Info<< nl
334  << "Initial feature set:" << nl
335  << " feature points : " << newSet.featurePoints().size() << nl
336  << " feature edges : " << newSet.featureEdges().size() << nl
337  << " of which" << nl
338  << " region edges : " << newSet.nRegionEdges() << nl
339  << " external edges : " << newSet.nExternalEdges() << nl
340  << " internal edges : " << newSet.nInternalEdges() << nl
341  << endl;
342 
343  boolList surfBaffleRegions(surf.patches().size(), false);
344 
345  wordList surfBaffleNames;
346  surfaceDict.readIfPresent("baffles", surfBaffleNames);
347 
348  forAll(surf.patches(), pI)
349  {
350  const word& name = surf.patches()[pI].name();
351 
352  if (findIndex(surfBaffleNames, name) != -1)
353  {
354  Info<< "Adding baffle region " << name << endl;
355  surfBaffleRegions[pI] = true;
356  }
357  }
358 
359  // Extracting and writing a extendedFeatureEdgeMesh
361  (
362  newSet,
363  runTime,
364  sFeatFileName + ".extendedFeatureEdgeMesh",
365  surfBaffleRegions
366  );
367 
368 
369  if (surfaceDict.isDict("addFeatures"))
370  {
371  const word addFeName = surfaceDict.subDict("addFeatures")["name"];
372  Info<< "Adding (without merging) features from " << addFeName
373  << nl << endl;
374 
375  extendedFeatureEdgeMesh addFeMesh
376  (
377  IOobject
378  (
379  addFeName,
380  runTime.time().constant(),
381  "extendedFeatureEdgeMesh",
382  runTime.time(),
385  )
386  );
387  Info<< "Read " << addFeMesh.name() << nl;
388  writeStats(addFeMesh, Info);
389 
390  feMesh.add(addFeMesh);
391  }
392 
393 
394  Info<< nl
395  << "Final feature set:" << nl;
396  writeStats(feMesh, Info);
397 
398  Info<< nl << "Writing extendedFeatureEdgeMesh to "
399  << feMesh.objectPath() << endl;
400 
401  mkDir(feMesh.path());
402 
403  if (writeObj)
404  {
405  feMesh.writeObj(feMesh.path()/surfFileName.lessExt().name());
406  }
407 
408  feMesh.write();
409 
410  // Write a featureEdgeMesh for backwards compatibility
411  featureEdgeMesh bfeMesh
412  (
413  IOobject
414  (
415  surfFileName.lessExt().name() + ".eMesh", // name
416  runTime.constant(), // instance
417  "triSurface",
418  runTime, // registry
421  false
422  ),
423  feMesh.points(),
424  feMesh.edges()
425  );
426 
427  Info<< nl << "Writing featureEdgeMesh to "
428  << bfeMesh.objectPath() << endl;
429 
430  bfeMesh.regIOobject::write();
431 
432 
433  // Find distance between close features
434  if (closeness)
435  {
436  Info<< nl << "Extracting internal and external closeness of "
437  << "surface." << endl;
438 
439  // Searchable triSurface
440  const triSurfaceMesh searchSurf
441  (
442  IOobject
443  (
444  sFeatFileName + ".closeness",
445  runTime.constant(),
446  "triSurface",
447  runTime
448  ),
449  surf
450  );
451 
452  {
453  Pair<tmp<triSurfaceScalarField>> closenessFields
454  (
455  searchSurf.extractCloseness()
456  );
457 
458  closenessFields.first()->write();
459  closenessFields.second()->write();
460 
461  if (writeVTK)
462  {
463  const faceList faces(searchSurf.faces());
464 
466  (
467  runTime.constantPath()/"triSurface",// outputDir
468  searchSurf.objectRegistry::name(), // surfaceName
469  searchSurf.points(),
470  faces,
471  "internalCloseness", // fieldName
472  closenessFields.first(),
473  false, // isNodeValues
474  true // verbose
475  );
476 
478  (
479  runTime.constantPath()/"triSurface",// outputDir
480  searchSurf.objectRegistry::name(), // surfaceName
481  searchSurf.points(),
482  faces,
483  "externalCloseness", // fieldName
484  closenessFields.second(),
485  false, // isNodeValues
486  true // verbose
487  );
488  }
489  }
490 
491  {
493  (
494  searchSurf.extractPointCloseness()
495  );
496 
497  closenessFields.first()->write();
498  closenessFields.second()->write();
499 
500  if (writeVTK)
501  {
502  const faceList faces(searchSurf.faces());
503  const Map<label>& meshPointMap = searchSurf.meshPointMap();
504 
506  internalClosenessPointField = closenessFields.first();
507 
509  externalClosenessPointField = closenessFields.second();
510 
511  scalarField internalCloseness(searchSurf.nPoints(), great);
512  scalarField externalCloseness(searchSurf.nPoints(), great);
513 
514  forAll(meshPointMap, pi)
515  {
516  internalCloseness[pi] =
517  internalClosenessPointField[meshPointMap[pi]];
518 
519  externalCloseness[pi] =
520  externalClosenessPointField[meshPointMap[pi]];
521  }
522 
524  (
525  runTime.constantPath()/"triSurface",// outputDir
526  searchSurf.objectRegistry::name(), // surfaceName
527  searchSurf.points(),
528  faces,
529  "internalPointCloseness", // fieldName
530  internalCloseness,
531  true, // isNodeValues
532  true // verbose
533  );
534 
536  (
537  runTime.constantPath()/"triSurface",// outputDir
538  searchSurf.objectRegistry::name(), // surfaceName
539  searchSurf.points(),
540  faces,
541  "externalPointCloseness", // fieldName
542  externalCloseness,
543  true, // isNodeValues
544  true // verbose
545  );
546  }
547  }
548  }
549 
550 
551  if (curvature)
552  {
553  Info<< nl << "Extracting curvature of surface at the points."
554  << endl;
555 
557  (
558  IOobject
559  (
560  sFeatFileName + ".curvature",
561  runTime.constant(),
562  "triSurface",
563  runTime
564  ),
565  surf,
566  dimLength,
567  surf.curvature()
568  );
569 
570  k.write();
571 
572  if (writeVTK)
573  {
575  (
576  runTime.constantPath()/"triSurface",// outputDir
577  sFeatFileName, // surfaceName
578  surf.points(),
579  faces,
580  "curvature", // fieldName
581  k,
582  true, // isNodeValues
583  true // verbose
584  );
585  }
586  }
587 
588 
589  if (featureProximity)
590  {
591  Info<< nl << "Extracting proximity of close feature points and "
592  << "edges to the surface" << endl;
593 
594  const scalar searchDistance =
595  readScalar(surfaceDict.lookup("maxFeatureProximity"));
596 
597  scalarField featureProximity(surf.size(), searchDistance);
598 
599  forAll(surf, fi)
600  {
601  const triPointRef& tri = surf[fi].tri(surf.points());
602  const point& triCentre = tri.circumCentre();
603 
604  const scalar radiusSqr = min
605  (
606  sqr(4*tri.circumRadius()),
607  sqr(searchDistance)
608  );
609 
610  pointIndexHitList hitList;
611 
612  feMesh.allNearestFeatureEdges(triCentre, radiusSqr, hitList);
613  featureProximity[fi] = min
614  (
615  feMesh.minDisconnectedDist(hitList),
616  featureProximity[fi]
617  );
618 
619  feMesh.allNearestFeaturePoints(triCentre, radiusSqr, hitList);
620  featureProximity[fi] = min
621  (
622  minDist(hitList),
623  featureProximity[fi]
624  );
625  }
626 
627  triSurfaceScalarField featureProximityField
628  (
629  IOobject
630  (
631  sFeatFileName + ".featureProximity",
632  runTime.constant(),
633  "triSurface",
634  runTime,
637  ),
638  surf,
639  dimLength,
640  featureProximity
641  );
642 
643  featureProximityField.write();
644 
645  if (writeVTK)
646  {
648  (
649  runTime.constantPath()/"triSurface",// outputDir
650  sFeatFileName, // surfaceName
651  surf.points(),
652  faces,
653  "featureProximity", // fieldName
654  featureProximity,
655  false, // isNodeValues
656  true // verbose
657  );
658  }
659  }
660 
661  Info<< endl;
662  }
663 
664  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
665  << " ClockTime = " << runTime.elapsedClockTime() << " s"
666  << nl << endl;
667 
668  Info<< "End\n" << endl;
669 
670  return 0;
671 }
672 
673 
674 // ************************************************************************* //
label nonFeatureStart() const
Return the index of the start of the non-feature points.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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
label multipleStart() const
Return the index of the start of the multiply-connected feature.
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
label openStart() const
Return the index of the start of the open feature edges.
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
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.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
label flatStart() const
Return the index of the start of the flat feature edges.
label internalStart() const
Return the index of the start of the internal feature edges.
A class for handling words, derived from string.
Definition: word.H:59
const edgeList & edges() const
Return edges.
Definition: edgeMeshI.H:68
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 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
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label mixedStart() const
Return the index of the start of the mixed type feature points.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
virtual void write(const fileName &outputDir, const fileName &surfaceName, const pointField &points, const faceList &faces, const bool verbose=false) const
Write single surface geometry to file.
const pointField & points() const
Return points.
Definition: edgeMeshI.H:62
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)
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
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
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.
label concaveStart() const
Return the index of the start of the concave feature points.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583