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