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-2021 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 "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 = dict.lookup<scalar>("includedAngle");
201 
203  (
204  surfaceFeatureSet
205  (
206  surfaceFileName,
207  surf,
208  dict,
209  includedAngle
210  )
211  );
212 
213  // Trim set
214  // ~~~~~~~~
215 
216  if (dict.isDict("trimFeatures"))
217  {
218  dictionary trimDict = dict.subDict("trimFeatures");
219 
220  scalar minLen =
221  trimDict.lookupOrAddDefault<scalar>("minLen", -great);
222 
223  label minElem = trimDict.lookupOrAddDefault<label>("minElem", 0);
224 
225  // Trim away small groups of features
226  if (minElem > 0 || minLen > 0)
227  {
228  Info<< "Removing features of length < "
229  << minLen << endl;
230  Info<< "Removing features with number of edges < "
231  << minElem << endl;
232 
233  set().trimFeatures(minLen, minElem, includedAngle);
234  }
235  }
236 
237 
238  // Subset
239  // ~~~~~~
240 
241  // Convert to marked edges, points
242  List<surfaceFeatures::edgeStatus> edgeStat(set().toStatus());
243 
244  if (dict.isDict("subsetFeatures"))
245  {
246  const dictionary& subsetDict = dict.subDict
247  (
248  "subsetFeatures"
249  );
250 
251  if (subsetDict.found("insideBox"))
252  {
253  treeBoundBox bb(subsetDict.lookup("insideBox")());
254 
255  Info<< "Selecting edges inside bb " << bb;
256  if (writeObj)
257  {
258  Info << " see insideBox.obj";
259  bb.writeOBJ("insideBox.obj");
260  }
261  Info<< endl;
262 
263  selectBox(surf, bb, true, edgeStat);
264  }
265  else if (subsetDict.found("outsideBox"))
266  {
267  treeBoundBox bb(subsetDict.lookup("outsideBox")());
268 
269  Info<< "Removing all edges inside bb " << bb;
270  if (writeObj)
271  {
272  Info<< " see outsideBox.obj" << endl;
273  bb.writeOBJ("outsideBox.obj");
274  }
275  Info<< endl;
276 
277  selectBox(surf, bb, false, edgeStat);
278  }
279 
280  const Switch nonManifoldEdges =
281  subsetDict.lookupOrDefault<Switch>("nonManifoldEdges", "yes");
282 
283  if (!nonManifoldEdges)
284  {
285  Info<< "Removing all non-manifold edges"
286  << " (edges with > 2 connected faces) unless they"
287  << " cross multiple regions" << endl;
288 
289  selectManifoldEdges(surf, 1e-5, includedAngle, edgeStat);
290  }
291 
292  const Switch openEdges =
293  subsetDict.lookupOrDefault<Switch>("openEdges", "yes");
294 
295  if (!openEdges)
296  {
297  Info<< "Removing all open edges"
298  << " (edges with 1 connected face)" << endl;
299 
300  forAll(edgeStat, edgei)
301  {
302  if (surf.edgeFaces()[edgei].size() == 1)
303  {
304  edgeStat[edgei] = surfaceFeatures::NONE;
305  }
306  }
307  }
308 
309  if (subsetDict.found("plane"))
310  {
311  const plane cutPlane(subsetDict.subDict("plane"));
312 
313  selectCutEdges(surf, cutPlane, edgeStat);
314 
315  Info<< "Only edges that intersect the plane with normal "
316  << cutPlane.normal()
317  << " and base point " << cutPlane.refPoint()
318  << " will be included as feature edges."<< endl;
319  }
320  }
321 
322 
323  surfaceFeatures newSet(surf);
324  newSet.setFromStatus(edgeStat, includedAngle);
325 
326  Info<< nl
327  << "Initial feature set:" << nl
328  << " feature points : " << newSet.featurePoints().size() << nl
329  << " feature edges : " << newSet.featureEdges().size() << nl
330  << " of which" << nl
331  << " region edges : " << newSet.nRegionEdges() << nl
332  << " external edges : " << newSet.nExternalEdges() << nl
333  << " internal edges : " << newSet.nInternalEdges() << nl
334  << endl;
335 
336  boolList surfBaffleRegions(surf.patches().size(), false);
337 
338  wordList surfBaffleNames;
339  dict.readIfPresent("baffles", surfBaffleNames);
340 
341  forAll(surf.patches(), pI)
342  {
343  const word& name = surf.patches()[pI].name();
344 
345  if (findIndex(surfBaffleNames, name) != -1)
346  {
347  Info<< "Adding baffle region " << name << endl;
348  surfBaffleRegions[pI] = true;
349  }
350  }
351 
352  // Extracting and writing a extendedFeatureEdgeMesh
354  (
355  newSet,
356  runTime,
357  sFeatFileName + ".extendedFeatureEdgeMesh",
358  surfBaffleRegions
359  );
360 
361 
362  if (dict.isDict("addFeatures"))
363  {
364  const word addFeName = dict.subDict("addFeatures")["name"];
365  Info<< "Adding (without merging) features from " << addFeName
366  << nl << endl;
367 
368  extendedFeatureEdgeMesh addFeMesh
369  (
370  IOobject
371  (
372  addFeName,
373  runTime.time().constant(),
374  "extendedFeatureEdgeMesh",
375  runTime.time(),
378  )
379  );
380  Info<< "Read " << addFeMesh.name() << nl;
381  addFeMesh.writeStats(Info);
382 
383  feMesh.add(addFeMesh);
384  }
385 
386 
387  Info<< nl
388  << "Final feature set:" << nl;
389  feMesh.writeStats(Info);
390 
391  Info<< nl << "Writing extendedFeatureEdgeMesh to "
392  << feMesh.localObjectPath() << endl;
393 
394  mkDir(feMesh.path());
395 
396  if (writeObj)
397  {
398  feMesh.writeObj
399  (
400  feMesh.path()/surfaceFileName.lessExt().name(),
401  verboseObj
402  );
403  }
404 
405  feMesh.write();
406 
407  // Write a featureEdgeMesh for backwards compatibility
408  featureEdgeMesh bfeMesh
409  (
410  IOobject
411  (
412  surfaceFileName.lessExt().name() + ".eMesh",
413  runTime.constant(),
415  runTime,
418  false
419  ),
420  feMesh.points(),
421  feMesh.edges()
422  );
423 
424  Info<< nl << "Writing featureEdgeMesh to "
425  << bfeMesh.localObjectPath() << endl;
426 
427  bfeMesh.regIOobject::write();
428 
429 
430  // Find distance between close features
431  if (dict.isDict("closeness"))
432  {
433  Info<< nl << "Extracting internal and external closeness of "
434  << "surface." << endl;
435 
436  const dictionary& closenessDict = dict.subDict("closeness");
437 
438  const Switch faceCloseness =
439  closenessDict.lookupOrDefault<Switch>("faceCloseness", "off");
440  const Switch pointCloseness =
441  closenessDict.lookupOrDefault<Switch>("pointCloseness", "off");
442 
443  const scalar internalAngleTolerance
444  (
445  closenessDict.lookupOrDefault<scalar>
446  (
447  "internalAngleTolerance", 80
448  )
449  );
450 
451  const scalar externalAngleTolerance
452  (
453  closenessDict.lookupOrDefault<scalar>
454  (
455  "externalAngleTolerance", 80
456  )
457  );
458 
459  // Searchable triSurface
460  const triSurfaceMesh searchSurf
461  (
462  IOobject
463  (
464  sFeatFileName + ".closeness",
465  runTime.constant(),
467  runTime
468  ),
469  surf
470  );
471 
472  if (faceCloseness)
473  {
474  Pair<tmp<triSurfaceScalarField>> closenessFields
475  (
476  searchSurf.extractCloseness
477  (
478  internalAngleTolerance,
479  externalAngleTolerance
480  )
481  );
482 
483  Info<< " writing "
484  << closenessFields.first()->name() << endl;
485  closenessFields.first()->write();
486 
487  Info<< " writing "
488  << closenessFields.second()->name() << endl;
489  closenessFields.second()->write();
490 
491  if (writeVTK)
492  {
493  const faceList faces(searchSurf.faces());
494 
496  (
497  runTime.path()
498  /runTime.constant()
500  searchSurf.objectRegistry::name(),
501  searchSurf.points(),
502  faces,
503  "internalCloseness", // fieldName
504  closenessFields.first(),
505  false // isNodeValues
506  );
507 
509  (
510  runTime.path()
511  /runTime.constant()
513  searchSurf.objectRegistry::name(),
514  searchSurf.points(),
515  faces,
516  "externalCloseness", // fieldName
517  closenessFields.second(),
518  false // isNodeValues
519  );
520  }
521  }
522 
523  if (pointCloseness)
524  {
526  (
527  searchSurf.extractPointCloseness
528  (
529  internalAngleTolerance,
530  externalAngleTolerance
531  )
532  );
533 
534  Info<< " writing "
535  << closenessFields.first()->name() << endl;
536  closenessFields.first()->write();
537 
538  Info<< " writing "
539  << closenessFields.second()->name() << endl;
540  closenessFields.second()->write();
541 
542  if (writeVTK)
543  {
544  const faceList faces(searchSurf.faces());
545  const Map<label>& meshPointMap = searchSurf.meshPointMap();
546 
548  internalClosenessPointField = closenessFields.first();
549 
551  externalClosenessPointField = closenessFields.second();
552 
553  scalarField internalCloseness(searchSurf.nPoints(), great);
554  scalarField externalCloseness(searchSurf.nPoints(), great);
555 
556  forAll(meshPointMap, pi)
557  {
558  internalCloseness[pi] =
559  internalClosenessPointField[meshPointMap[pi]];
560 
561  externalCloseness[pi] =
562  externalClosenessPointField[meshPointMap[pi]];
563  }
564 
566  (
567  runTime.path()
568  /runTime.constant()
570  searchSurf.objectRegistry::name(),
571  searchSurf.points(),
572  faces,
573  "internalPointCloseness", // fieldName
574  internalCloseness,
575  true // isNodeValues
576  );
577 
579  (
580  runTime.path()
581  /runTime.constant()
583  searchSurf.objectRegistry::name(),
584  searchSurf.points(),
585  faces,
586  "externalPointCloseness", // fieldName
587  externalCloseness,
588  true // isNodeValues
589  );
590  }
591  }
592  }
593 
594 
595  if (curvature)
596  {
597  Info<< nl << "Extracting curvature of surface at the points."
598  << endl;
599 
601  (
602  IOobject
603  (
604  sFeatFileName + ".curvature",
605  runTime.constant(),
607  runTime
608  ),
609  surf,
610  dimLength,
611  surf.curvature()
612  );
613 
614  k.write();
615 
616  if (writeVTK)
617  {
619  (
620  runTime.path()
621  /runTime.constant()
623  sFeatFileName,
624  surf.points(),
625  faces,
626  "curvature", // fieldName
627  k,
628  true // isNodeValues
629  );
630  }
631  }
632 
633 
634  if (featureProximity)
635  {
636  Info<< nl << "Extracting proximity of close feature points and "
637  << "edges to the surface" << endl;
638 
639  const scalar searchDistance =
640  dict.lookup<scalar>("maxFeatureProximity");
641 
642  scalarField featureProximity(surf.size(), searchDistance);
643 
644  forAll(surf, fi)
645  {
646  const triPointRef& tri = surf[fi].tri(surf.points());
647  const point& triCentre = tri.circumCentre();
648 
649  const scalar radiusSqr = min
650  (
651  sqr(4*tri.circumRadius()),
652  sqr(searchDistance)
653  );
654 
655  pointIndexHitList hitList;
656 
657  feMesh.allNearestFeatureEdges(triCentre, radiusSqr, hitList);
658  featureProximity[fi] = min
659  (
660  feMesh.minDisconnectedDist(hitList),
661  featureProximity[fi]
662  );
663 
664  feMesh.allNearestFeaturePoints(triCentre, radiusSqr, hitList);
665  featureProximity[fi] = min
666  (
667  minDist(hitList),
668  featureProximity[fi]
669  );
670  }
671 
672  triSurfaceScalarField featureProximityField
673  (
674  IOobject
675  (
676  sFeatFileName + ".featureProximity",
677  runTime.constant(),
679  runTime,
682  ),
683  surf,
684  dimLength,
685  featureProximity
686  );
687 
688  featureProximityField.write();
689 
690  if (writeVTK)
691  {
693  (
694  runTime.path()
695  /runTime.constant()
697  sFeatFileName,
698  surf.points(),
699  faces,
700  "featureProximity", // fieldName
701  featureProximity,
702  false // isNodeValues
703  );
704  }
705  }
706 
707  Info<< endl;
708  }
709 
710 
711  void extractFeatures
712  (
713  const fileNameList& surfaceFileNames,
714  const Time& runTime,
715  const dictionary& dict
716  )
717  {
718  forAll(surfaceFileNames, i)
719  {
720  extractFeatures(surfaceFileNames[i], runTime, dict);
721  }
722  }
723 }
724 
725 
726 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
727 
728 int main(int argc, char *argv[])
729 {
731  (
732  "extract and write surface features to file"
733  );
735 
736  #include "addDictOption.H"
737 
738  #include "setRootCase.H"
739  #include "createTime.H"
740 
741  const dictionary dict(systemDict("surfaceFeaturesDict", args, runTime));
742 
743  if (dict.found("surfaces"))
744  {
745  extractFeatures
746  (
747  fileNameList(dict.lookup("surfaces")),
748  runTime,
749  dict
750  );
751  }
752  else
753  {
754  forAllConstIter(dictionary, dict, iter)
755  {
756  if (!iter().isDict())
757  {
758  continue;
759  }
760 
761  extractFeatures
762  (
763  fileNameList(iter().dict().lookup("surfaces")),
764  runTime,
765  iter().dict()
766  );
767  }
768  }
769 
770  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
771  << " ClockTime = " << runTime.elapsedClockTime() << " s"
772  << nl << endl;
773 
774  Info<< "End\n" << endl;
775 
776  return 0;
777 }
778 
779 
780 
781 // ************************************************************************* //
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#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
virtual Ostream & write(const char)=0
Write character.
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:156
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:164
tmp< scalarField > curvature() const
Return the curvature of surface at the points.
Definition: triSurface.C:1144
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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:175
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
T & first()
Return the first element of the list.
Definition: UListI.H:114
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
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:936
const dimensionSet dimLength
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
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:1353
static const word & geometryDir()
Return the geometry directory name.
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
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.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
const labelListList & edgeFaces() const
Return edge-face addressing.
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:260
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:285
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
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 with support for writing ASCII or binary.
fileName lessExt() const
Return file name without extension (part before last .)
Definition: fileName.C:284
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:1131
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
Triangulated surface description with patch information.
Definition: triSurface.H:66
Foam::argList args(argc, argv)
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
fileName path() const
Return path.
Definition: TimePaths.H:139
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
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:844