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 "vtkWritePolyData.H"
38 #include "systemDict.H"
39 
40 using namespace Foam;
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  autoPtr<surfaceFeatures> extractFromFile
47  (
48  const fileName& featureEdgeFile,
49  const triSurface& surf,
50  const Switch& geometricTestOnly
51  )
52  {
53  edgeMesh eMesh(featureEdgeFile);
54 
55  // Sometimes duplicate edges are present. Remove them.
56  eMesh.mergeEdges();
57 
58  Info<< nl << "Reading existing feature edges from file "
59  << featureEdgeFile << nl
60  << "Selecting edges purely based on geometric tests: "
61  << geometricTestOnly.asText() << endl;
62 
64  (
65  new surfaceFeatures
66  (
67  surf,
68  eMesh.points(),
69  eMesh.edges(),
70  1e-6,
71  geometricTestOnly
72  )
73  );
74  }
75 
76 
77  autoPtr<surfaceFeatures> extractFromSurface
78  (
79  const triSurface& surf,
80  const Switch& geometricTestOnly,
81  const scalar includedAngle
82  )
83  {
84  Info<< nl
85  << "Constructing feature set from included angle "
86  << includedAngle << nl
87  << "Selecting edges purely based on geometric tests: "
88  << geometricTestOnly.asText() << endl;
89 
91  (
92  new surfaceFeatures
93  (
94  surf,
95  includedAngle,
96  0,
97  0,
98  geometricTestOnly
99  )
100  );
101  }
102 
103 
104  autoPtr<surfaceFeatures> surfaceFeatureSet
105  (
106  const fileName& surfaceFileName,
107  const triSurface& surf,
108  const dictionary& dict,
109  const scalar includedAngle
110  )
111  {
112  const Switch geometricTestOnly = dict.lookupOrDefault<Switch>
113  (
114  "geometricTestOnly",
115  "no"
116  );
117 
118  if (dict.found("files"))
119  {
120  HashTable<fileName, fileName> fileNames(dict.lookup("files"));
121 
122  if (fileNames.found(surfaceFileName))
123  {
124  return extractFromFile
125  (
126  fileNames[surfaceFileName],
127  surf,
128  geometricTestOnly
129  );
130  }
131  else
132  {
133  return extractFromSurface
134  (
135  surf,
136  geometricTestOnly,
137  includedAngle
138  );
139  }
140  }
141  else
142  {
143  return extractFromSurface
144  (
145  surf,
146  geometricTestOnly,
147  includedAngle
148  );
149  }
150  }
151 
152 
153  void extractFeatures
154  (
155  const fileName& surfaceFileName,
156  const Time& runTime,
157  const dictionary& dict
158  )
159  {
160  const fileName sFeatFileName = surfaceFileName.lessExt().name();
161 
162  Info<< "Surface : " << surfaceFileName << nl << endl;
163 
164  const Switch writeVTK =
165  dict.lookupOrDefault<Switch>("writeVTK", "off");
166  const Switch writeObj =
167  dict.lookupOrDefault<Switch>("writeObj", "off");
168  const Switch verboseObj =
169  dict.lookupOrDefault<Switch>("verboseObj", "off");
170 
171  const Switch curvature =
172  dict.lookupOrDefault<Switch>("curvature", "off");
173  const Switch featureProximity =
174  dict.lookupOrDefault<Switch>("featureProximity", "off");
175 
176 
177  Info<< nl << "Feature line extraction is only valid on closed manifold "
178  << "surfaces." << endl;
179 
180  // Read
181  // ~~~~
182 
183  triSurface surf
184  (
185  runTime.path()
186  /runTime.constant()
188  /surfaceFileName
189  );
190 
191  Info<< "Statistics:" << endl;
192  surf.writeStats(Info);
193  Info<< endl;
194 
195  const faceList faces(surf.faces());
196 
197  // Either construct features from surface & featureAngle or read set.
198  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
199 
200  const scalar includedAngle = 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.relativeObjectPath() << 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.relativeObjectPath() << endl;
426 
427  bfeMesh.regIOobject::write();
428 
429  // Data to write out in VTK format
430  wordList writeVTKFieldNames;
431  boolList writeVTKFieldIsPointValues;
432  #define DeclareWriteVTKFieldTypeValues(Type, nullArg) \
433  PtrList<const Field<Type>> writeVTKField##Type##Values;
434  FOR_ALL_FIELD_TYPES(DeclareWriteVTKFieldTypeValues);
435  #undef DeclareWriteVTKFieldTypeValues
436 
437  // Find distance between close features
438  if (dict.isDict("closeness"))
439  {
440  Info<< nl << "Extracting internal and external closeness of "
441  << "surface." << endl;
442 
443  const dictionary& closenessDict = dict.subDict("closeness");
444 
445  const Switch faceCloseness =
446  closenessDict.lookupOrDefault<Switch>("faceCloseness", "off");
447  const Switch pointCloseness =
448  closenessDict.lookupOrDefault<Switch>("pointCloseness", "off");
449 
450  const scalar internalAngleTolerance
451  (
452  closenessDict.lookupOrDefault<scalar>
453  (
454  "internalAngleTolerance", 80
455  )
456  );
457 
458  const scalar externalAngleTolerance
459  (
460  closenessDict.lookupOrDefault<scalar>
461  (
462  "externalAngleTolerance", 80
463  )
464  );
465 
466  // Searchable triSurface
467  const triSurfaceMesh searchSurf
468  (
469  IOobject
470  (
471  sFeatFileName + ".closeness",
472  runTime.constant(),
474  runTime
475  ),
476  surf
477  );
478 
479  if (faceCloseness)
480  {
481  Pair<tmp<triSurfaceScalarField>> closenessFields
482  (
483  searchSurf.extractCloseness
484  (
485  internalAngleTolerance,
486  externalAngleTolerance
487  )
488  );
489 
490  Info<< " writing "
491  << closenessFields.first()->name() << endl;
492  closenessFields.first()->write();
493 
494  Info<< " writing "
495  << closenessFields.second()->name() << endl;
496  closenessFields.second()->write();
497 
498  if (writeVTK)
499  {
500  writeVTKFieldNames.append("internalCloseness");
501  writeVTKFieldIsPointValues.append(false);
502  writeVTKFieldscalarValues.append
503  (
504  new scalarField(closenessFields.first())
505  );
506 
507  writeVTKFieldNames.append("externalCloseness");
508  writeVTKFieldIsPointValues.append(false);
509  writeVTKFieldscalarValues.append
510  (
511  new scalarField(closenessFields.second())
512  );
513  }
514  }
515 
516  if (pointCloseness)
517  {
519  (
520  searchSurf.extractPointCloseness
521  (
522  internalAngleTolerance,
523  externalAngleTolerance
524  )
525  );
526 
527  Info<< " writing "
528  << closenessFields.first()->name() << endl;
529  closenessFields.first()->write();
530 
531  Info<< " writing "
532  << closenessFields.second()->name() << endl;
533  closenessFields.second()->write();
534 
535  if (writeVTK)
536  {
537  const faceList faces(searchSurf.faces());
538  const Map<label>& meshPointMap = searchSurf.meshPointMap();
539 
541  internalClosenessPointField = closenessFields.first();
542 
544  externalClosenessPointField = closenessFields.second();
545 
546  scalarField internalCloseness(searchSurf.nPoints(), great);
547  scalarField externalCloseness(searchSurf.nPoints(), great);
548 
549  forAll(meshPointMap, pi)
550  {
551  internalCloseness[pi] =
552  internalClosenessPointField[meshPointMap[pi]];
553 
554  externalCloseness[pi] =
555  externalClosenessPointField[meshPointMap[pi]];
556  }
557 
558  writeVTKFieldNames.append("internalPointCloseness");
559  writeVTKFieldIsPointValues.append(true);
560  writeVTKFieldscalarValues.append
561  (
562  new scalarField(internalCloseness)
563  );
564 
565  writeVTKFieldNames.append("externalPointCloseness");
566  writeVTKFieldIsPointValues.append(true);
567  writeVTKFieldscalarValues.append
568  (
569  new scalarField(externalCloseness)
570  );
571  }
572  }
573  }
574 
575 
576  if (curvature)
577  {
578  Info<< nl << "Extracting curvature of surface at the points."
579  << endl;
580 
582  (
583  IOobject
584  (
585  sFeatFileName + ".curvature",
586  runTime.constant(),
588  runTime
589  ),
590  surf,
591  dimLength,
592  surf.curvature()
593  );
594 
595  k.write();
596 
597  if (writeVTK)
598  {
599  writeVTKFieldNames.append("curvature");
600  writeVTKFieldIsPointValues.append(true);
601  writeVTKFieldscalarValues.append(new scalarField(k));
602  }
603  }
604 
605 
606  if (featureProximity)
607  {
608  Info<< nl << "Extracting proximity of close feature points and "
609  << "edges to the surface" << endl;
610 
611  const scalar searchDistance =
612  dict.lookup<scalar>("maxFeatureProximity");
613 
614  scalarField featureProximity(surf.size(), searchDistance);
615 
616  forAll(surf, fi)
617  {
618  const triPointRef& tri = surf[fi].tri(surf.points());
619  const point& triCentre = tri.circumCentre();
620 
621  const scalar radiusSqr = min
622  (
623  sqr(4*tri.circumRadius()),
624  sqr(searchDistance)
625  );
626 
627  pointIndexHitList hitList;
628 
629  feMesh.allNearestFeatureEdges(triCentre, radiusSqr, hitList);
630  featureProximity[fi] = min
631  (
632  feMesh.minDisconnectedDist(hitList),
633  featureProximity[fi]
634  );
635 
636  feMesh.allNearestFeaturePoints(triCentre, radiusSqr, hitList);
637  featureProximity[fi] = min
638  (
639  minDist(hitList),
640  featureProximity[fi]
641  );
642  }
643 
644  triSurfaceScalarField featureProximityField
645  (
646  IOobject
647  (
648  sFeatFileName + ".featureProximity",
649  runTime.constant(),
651  runTime,
654  ),
655  surf,
656  dimLength,
657  featureProximity
658  );
659 
660  featureProximityField.write();
661 
662  if (writeVTK)
663  {
664  writeVTKFieldNames.append("featureProximity");
665  writeVTKFieldIsPointValues.append(false);
666  writeVTKFieldscalarValues.append
667  (
668  new scalarField(featureProximity)
669  );
670  }
671  }
672 
673  if (writeVTK)
674  {
675  #define WriteVTKResizeFieldTypeValues(Type, nullArg) \
676  writeVTKField##Type##Values.resize(writeVTKFieldNames.size());
677  FOR_ALL_FIELD_TYPES(WriteVTKResizeFieldTypeValues)
678  #undef WriteVTKResizeFieldTypeValues
679 
681  (
682  runTime.path()
683  /runTime.constant()
685  /sFeatFileName + "Features.vtk",
686  sFeatFileName,
687  runTime.writeFormat() == IOstream::BINARY,
688  surf.points(),
689  labelList(),
690  labelListList(),
691  faces,
692  writeVTKFieldNames,
693  writeVTKFieldIsPointValues,
694  UPtrList<const Field<label>>(writeVTKFieldNames.size())
695  #define WriteVTKFieldTypeValuesParameter(Type, nullArg) \
696  , UPtrList<const Field<Type>>(writeVTKField##Type##Values)
697  FOR_ALL_FIELD_TYPES(WriteVTKFieldTypeValuesParameter)
698  #undef WriteVTKFieldTypeValuesParameter
699  );
700  }
701 
702  Info<< endl;
703  }
704 
705 
706  void extractFeatures
707  (
708  const fileNameList& surfaceFileNames,
709  const Time& runTime,
710  const dictionary& dict
711  )
712  {
713  forAll(surfaceFileNames, i)
714  {
715  extractFeatures(surfaceFileNames[i], runTime, dict);
716  }
717  }
718 }
719 
720 
721 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
722 
723 int main(int argc, char *argv[])
724 {
726  (
727  "extract and write surface features to file"
728  );
730 
731  #include "addDictOption.H"
732 
733  #include "setRootCase.H"
734  #include "createTime.H"
735 
736  const dictionary dict(systemDict("surfaceFeaturesDict", args, runTime));
737 
738  if (dict.found("surfaces"))
739  {
740  extractFeatures
741  (
742  fileNameList(dict.lookup("surfaces")),
743  runTime,
744  dict
745  );
746  }
747  else
748  {
749  forAllConstIter(dictionary, dict, iter)
750  {
751  if (!iter().isDict())
752  {
753  continue;
754  }
755 
756  extractFeatures
757  (
758  fileNameList(iter().dict().lookup("surfaces")),
759  runTime,
760  iter().dict()
761  );
762  }
763  }
764 
765  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
766  << " ClockTime = " << runTime.elapsedClockTime() << " s"
767  << nl << endl;
768 
769  Info<< "End\n" << endl;
770 
771  return 0;
772 }
773 
774 
775 
776 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#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.
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
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
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
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:69
scalar minDist(const List< pointIndexHit > &hitList)
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:956
const dimensionSet dimLength
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
A surface geometry formed of discrete facets, e.g. triangles and/or quadrilaterals, defined in a file using formats such as Wavefront OBJ, or stereolithography STL.
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
static const word & geometryDir()
Return the geometry directory name.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
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.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const labelListList & edgeFaces() const
Return edge-face addressing.
List< label > labelList
A List of labels.
Definition: labelList.H:56
scalar circumRadius() const
Return circum-radius.
Definition: triangleI.H:146
An STL-conforming hash table.
Definition: HashTable.H:61
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
static const char nl
Definition: Ostream.H:260
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:293
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
void write(const fileName &file, const word &title, const bool binary, const PointField &points, const VertexList &vertices, const LineList &lines, const FaceList &faces, const wordList &fieldNames, const boolList &fieldIsPointValues, const UPtrList< const Field< label >> &fieldLabelValues #define FieldTypeValuesConstArg(Type, nullArg))
Write VTK polygonal data to a file. Takes a PtrList of fields of labels and.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
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:53
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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
FOR_ALL_FIELD_TYPES(DefineFvWallInfoType)
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
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:98
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:864