surfaceFeatureExtract.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  surfaceFeatureExtract
29 
30 Group
31  grpSurfaceUtilities
32 
33 Description
34  Extracts and writes surface features to file. All but the basic feature
35  extraction is a work-in-progress.
36 
37  The extraction process is driven by the \a system/surfaceFeatureExtractDict
38  dictionary, but the \a -dict option can be used to define an alternative
39  location.
40 
41  The \a system/surfaceFeatureExtractDict dictionary contains entries
42  for each extraction process.
43  The name of the individual dictionary is used to load the input surface
44  (found under \a constant/triSurface) and also as the basename for the
45  output.
46 
47  If the \c surfaces entry is present in a sub-dictionary, it has absolute
48  precedence over a surface name deduced from the dictionary name.
49  If the dictionary name itself does not have an extension, the \c surfaces
50  entry becomes mandatory since in this case the dictionary name cannot
51  represent an input surface file (ie, there is no file extension).
52  The \c surfaces entry is a wordRe list, which allows loading and
53  combining of multiple surfaces. Any exactly specified surface names must
54  exist, but surfaces selected via regular expressions need not exist.
55  The selection mechanism preserves order and is without duplicates.
56  For example,
57  \verbatim
58  dictName
59  {
60  surfaces (surface1.stl "other.*" othersurf.obj);
61  ...
62  }
63  \endverbatim
64 
65  When loading surfaces, the points/faces/regions of each surface are
66  normally offset to create an aggregated surface. No merging of points
67  or faces is done. The optional entry \c loadingOption can be used to
68  adjust the treatment of the regions when loading single or multiple files,
69  with selections according to the Foam::triSurfaceLoader::loadingOption
70  enumeration.
71  \verbatim
72  dictName
73  {
74  // Optional treatment of surface regions when loading
75  // (single, file, offset, merge)
76  loadingOption file;
77  ...
78  }
79  \endverbatim
80  The \c loadingOption is primarily used in combination with the
81  \c intersectionMethod (specifically its \c region option).
82  The default \c loadingOption is normally \c offset,
83  but this changes to \c file if the \c intersectionMethod
84  \c region is being used.
85 
86  Once surfaces have been loaded, the first stage is to extract
87  features according to the specified \c extractionMethod with values
88  as per the following table:
89  \table
90  extractionMethod | Description
91  none | No feature extraction
92  extractFromFile | Load features from the file named in featureEdgeFile
93  extractFromSurface | Extract features from surface geometry
94  \endtable
95 
96  There are a few entries that influence the extraction behaviour:
97  \verbatim
98  // File to use for extractFromFile input
99  featureEdgeFile "FileName"
100 
101  // Mark edges whose adjacent surface normals are at an angle less
102  // than includedAngle as features
103  // - 0 : selects no edges
104  // - 180: selects all edges
105  includedAngle 120;
106 
107  // Do not mark region edges
108  geometricTestOnly yes;
109  \endverbatim
110 
111  This initial set of edges can be trimmed:
112  \verbatim
113  trimFeatures
114  {
115  // Remove features with fewer than the specified number of edges
116  minElem 0;
117 
118  // Remove features shorter than the specified cumulative length
119  minLen 0.0;
120  }
121  \endverbatim
122 
123  and subsetted
124  \verbatim
125  subsetFeatures
126  {
127  // Use a plane to select feature edges (normal)(basePoint)
128  // Only keep edges that intersect the plane
129  plane (1 0 0)(0 0 0);
130 
131  // Select feature edges using a box // (minPt)(maxPt)
132  // Only keep edges inside the box:
133  insideBox (0 0 0)(1 1 1);
134 
135  // Only keep edges outside the box:
136  outsideBox (0 0 0)(1 1 1);
137 
138  // Keep nonManifold edges (edges with >2 connected faces where
139  // the faces form more than two different normal planes)
140  nonManifoldEdges yes;
141 
142  // Keep open edges (edges with 1 connected face)
143  openEdges yes;
144  }
145  \endverbatim
146 
147  Subsequently, additional features can be added from another file:
148  \verbatim
149  addFeatures
150  {
151  // Add (without merging) another extendedFeatureEdgeMesh
152  name axZ.extendedFeatureEdgeMesh;
153  }
154  \endverbatim
155 
156  The intersectionMethod provides a final means of adding additional
157  features. These are loosely termed "self-intersection", since it
158  detects the face/face intersections of the loaded surface or surfaces.
159 
160  \table
161  intersectionMethod | Description
162  none | Do nothing
163  self | All face/face intersections
164  region | Limit face/face intersections to those between different regions.
165  \endtable
166  The optional \c tolerance tuning parameter is available for handling
167  the face/face intersections, but should normally not be touched.
168 
169  As well as the normal extendedFeatureEdgeMesh written,
170  other items can be selected with boolean switches:
171 
172  \table
173  Output option | Description
174  closeness | Output the closeness of surface elements to other surface elements.
175  curvature | Output surface curvature
176  featureProximity | Output the proximity of feature points and edges to another
177  writeObj | Write features to OBJ format for postprocessing
178  writeVTK | Write closeness/curvature/proximity fields as VTK for postprocessing
179  \endtable
180 
181 Note
182  Although surfaceFeatureExtract can do many things, there are still a fair
183  number of corner cases where it may not produce the desired result.
184 \*---------------------------------------------------------------------------*/
185 
186 #include "argList.H"
187 #include "Time.H"
188 #include "triSurface.H"
189 #include "triSurfaceTools.H"
190 #include "edgeMeshTools.H"
192 #include "surfaceIntersection.H"
193 #include "featureEdgeMesh.H"
194 #include "extendedFeatureEdgeMesh.H"
195 #include "treeBoundBox.H"
196 #include "meshTools.H"
197 #include "OBJstream.H"
198 #include "triSurfaceMesh.H"
199 #include "foamVtkSurfaceWriter.H"
200 #include "unitConversion.H"
201 #include "plane.H"
202 #include "point.H"
203 #include "triSurfaceLoader.H"
204 
205 using namespace Foam;
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 int main(int argc, char *argv[])
210 {
212  (
213  "Extract and write surface feature lines to file.\n"
214  "Feature line extraction only valid on closed manifold surfaces."
215  );
216 
218  argList::noFunctionObjects(); // Never use function objects
219 
221  (
222  "dict",
223  "file",
224  "Read surfaceFeatureExtractDict from specified location"
225  );
226 
227  #include "setRootCase.H"
228  #include "createTime.H"
229 
230  Info<< nl
231  << "Note: "
232  << "Feature line extraction only valid on closed manifold surfaces"
233  << nl << nl;
234 
235  const word dictName("surfaceFeatureExtractDict");
237 
238  Info<< "Reading " << dictIO.name() << nl << endl;
239  const IOdictionary dict(dictIO);
240 
241  // Loader for available triSurface surface files
242  triSurfaceLoader loader(runTime);
243 
244  // Where to write VTK output files
245  const fileName vtkOutputDir = runTime.constantPath()/"triSurface";
246 
247  for (const entry& dEntry : dict)
248  {
249  if (!dEntry.isDict()) // dictionary entries only
250  {
251  continue;
252  }
253 
254  // Note - do not check for shell meta-characters
255  // - user responsibility
256 
257  const word& dictName = dEntry.keyword();
258  const dictionary& surfaceDict = dEntry.dict();
259 
260  if (!surfaceDict.found("extractionMethod"))
261  {
262  // Insist on an extractionMethod
263  continue;
264  }
265 
266  // The output name based in dictionary name (without extensions)
267  const word outputName = dictName.lessExt();
268 
271  (
272  surfaceDict
273  );
274 
275  // We don't needs the intersectionMethod yet, but can use it
276  // for setting a reasonable loading option
277  const surfaceIntersection::intersectionType selfIntersect =
279  (
280  "intersectionMethod",
281  surfaceDict,
283  );
284 
285  const Switch writeObj("writeObj", surfaceDict, Switch::OFF);
286 
287  const Switch writeVTK("writeVTK", surfaceDict, Switch::OFF);
288 
289  // The "surfaces" entry is normally optional, but make it mandatory
290  // if the dictionary name doesn't have an extension
291  // (ie, probably not a surface filename at all).
292  // If it is missing, this will fail nicely with an appropriate error
293  // message.
294  if (surfaceDict.found("surfaces") || !dictName.hasExt())
295  {
296  loader.select(surfaceDict.get<wordRes>("surfaces"));
297  }
298  else
299  {
300  loader.select(dictName);
301  }
302 
303  // DebugVar(loader.available());
304  // DebugVar(outputName);
305 
306  if (loader.selected().empty())
307  {
309  << "No surfaces specified/found for entry: "
310  << dictName << exit(FatalError);
311  }
312 
313  Info<< "Surfaces : ";
314  if (loader.selected().size() == 1)
315  {
316  Info<< loader.selected().first() << nl;
317  }
318  else
319  {
320  Info<< flatOutput(loader.selected()) << nl;
321  }
322  Info<< "Output : " << outputName << nl;
323 
324  // Loading option - default depends on context
325  triSurfaceLoader::loadingOption loadingOption =
327  (
328  "loadingOption",
329  surfaceDict,
330  (
331  selfIntersect == surfaceIntersection::SELF_REGION
334  )
335  );
336 
337  Info<<"Load options : "
338  << triSurfaceLoader::loadingOptionNames[loadingOption] << nl
339  << "Write options:"
340  << " writeObj=" << writeObj
341  << " writeVTK=" << writeVTK << nl;
342 
343  scalar scaleFactor = -1;
344  // Allow rescaling of the surface points (eg, mm -> m)
345  if (surfaceDict.readIfPresent("scale", scaleFactor) && scaleFactor > 0)
346  {
347  Info<<"Scaling : " << scaleFactor << nl;
348  }
349 
350  // Load a single file, or load and combine multiple selected files
351  autoPtr<triSurface> surfPtr = loader.load(loadingOption, scaleFactor);
352  if (!surfPtr || surfPtr->empty())
353  {
355  << "Problem loading surface(s) for entry: "
356  << dictName << exit(FatalError);
357  }
358 
359  triSurface surf = *surfPtr;
360 
361  Info<< nl
362  << "Statistics:" << nl;
363  surf.writeStats(Info);
364 
365  // Need a copy as plain faces if outputting VTK format
366  faceList faces;
367  if (writeVTK)
368  {
369  faces.setSize(surf.size());
370  forAll(surf, fi)
371  {
372  faces[fi] = surf[fi];
373  }
374  }
375 
376  //
377  // Extract features using the preferred extraction method
378  //
379  autoPtr<surfaceFeatures> features = extractor().features(surf);
380 
381  // Trim set
382  // ~~~~~~~~
383 
384  // Option: "trimFeatures" (dictionary)
385  if (surfaceDict.isDict("trimFeatures"))
386  {
387  const dictionary& trimDict = surfaceDict.subDict("trimFeatures");
388 
389  const scalar minLen =
390  trimDict.getOrDefault<scalar>("minLen", 0);
391  const label minElem =
392  trimDict.getOrDefault<label>("minElem", 0);
393 
394  // Trim away small groups of features
395  if (minLen > 0 || minElem > 0)
396  {
397  if (minLen > 0)
398  {
399  Info<< "Removing features of length < "
400  << minLen << endl;
401  }
402  if (minElem > 0)
403  {
404  Info<< "Removing features with number of edges < "
405  << minElem << endl;
406  }
407 
408  features().trimFeatures
409  (
410  minLen, minElem, extractor().includedAngle()
411  );
412  }
413  }
414 
415  // Subset
416  // ~~~~~~
417 
418  // Convert to marked edges, points
419  List<surfaceFeatures::edgeStatus> edgeStat(features().toStatus());
420 
421  // Option: "subsetFeatures" (dictionary)
422  if (surfaceDict.isDict("subsetFeatures"))
423  {
424  const dictionary& subsetDict = surfaceDict.subDict
425  (
426  "subsetFeatures"
427  );
428 
429  // Suboption: "insideBox"
430  if (subsetDict.found("insideBox"))
431  {
432  treeBoundBox bb(subsetDict.lookup("insideBox")());
433 
434  Info<< "Subset edges inside box " << bb << endl;
435  features().subsetBox(edgeStat, bb);
436 
437  {
438  OBJstream os("subsetBox.obj");
439 
440  Info<< "Dumping bounding box " << bb
441  << " as lines to obj file "
442  << os.name() << endl;
443 
444  os.write(bb);
445  }
446  }
447  // Suboption: "outsideBox"
448  else if (subsetDict.found("outsideBox"))
449  {
450  treeBoundBox bb(subsetDict.lookup("outsideBox")());
451 
452  Info<< "Exclude edges outside box " << bb << endl;
453  features().excludeBox(edgeStat, bb);
454 
455  {
456  OBJstream os("deleteBox.obj");
457 
458  Info<< "Dumping bounding box " << bb
459  << " as lines to obj file "
460  << os.name() << endl;
461 
462  os.write(bb);
463  }
464  }
465 
466  // Suboption: "nonManifoldEdges" (false: remove non-manifold edges)
467  if (!subsetDict.getOrDefault("nonManifoldEdges", true))
468  {
469  Info<< "Removing all non-manifold edges"
470  << " (edges with > 2 connected faces) unless they"
471  << " cross multiple regions" << endl;
472 
473  features().checkFlatRegionEdge
474  (
475  edgeStat,
476  1e-5, // tol
477  extractor().includedAngle()
478  );
479  }
480 
481  // Suboption: "openEdges" (false: remove open edges)
482  if (!subsetDict.getOrDefault("openEdges", true))
483  {
484  Info<< "Removing all open edges"
485  << " (edges with 1 connected face)" << endl;
486 
487  features().excludeOpen(edgeStat);
488  }
489 
490  // Suboption: "plane"
491  if (subsetDict.found("plane"))
492  {
493  plane cutPlane(subsetDict.lookup("plane")());
494 
495  Info<< "Only include feature edges that intersect the plane"
496  << " with normal " << cutPlane.normal()
497  << " and origin " << cutPlane.origin() << endl;
498 
499  features().subsetPlane(edgeStat, cutPlane);
500  }
501  }
502 
503  surfaceFeatures newSet(surf);
504  newSet.setFromStatus(edgeStat, extractor().includedAngle());
505 
506  Info<< nl << "Initial ";
507  newSet.writeStats(Info);
508 
509  boolList surfBaffleRegions(surf.patches().size(), false);
510  if (surfaceDict.found("baffles"))
511  {
512  wordRes baffleSelect(surfaceDict.get<wordRes>("baffles"));
513 
514  wordList patchNames(surf.patches().size());
515  forAll(surf.patches(), patchi)
516  {
517  patchNames[patchi] = surf.patches()[patchi].name();
518  }
519 
520  labelList indices(baffleSelect.matching(patchNames));
521 
522  for (const label patchId : indices)
523  {
524  surfBaffleRegions[patchId] = true;
525  }
526 
527  if (indices.size())
528  {
529  Info<< "Adding " << indices.size() << " baffle regions: (";
530 
531  forAll(surfBaffleRegions, patchi)
532  {
533  if (surfBaffleRegions[patchi])
534  {
535  Info<< ' ' << patchNames[patchi];
536  }
537  }
538  Info<< " )" << nl << nl;
539  }
540  }
541 
542  // Extracting and writing a extendedFeatureEdgeMesh
544  (
545  newSet,
546  runTime,
547  outputName + ".extendedFeatureEdgeMesh",
548  surfBaffleRegions
549  );
550 
551 
552  if (surfaceDict.isDict("addFeatures"))
553  {
554  const word addFeName
555  (
556  surfaceDict.subDict("addFeatures").get<word>("name")
557  );
558 
559  Info<< "Adding (without merging) features from " << addFeName
560  << nl << endl;
561 
562  extendedFeatureEdgeMesh addFeMesh
563  (
564  IOobject
565  (
566  addFeName,
567  runTime.time().constant(),
568  "extendedFeatureEdgeMesh",
569  runTime.time(),
572  )
573  );
574  Info<< "Read " << addFeMesh.name() << nl;
575  edgeMeshTools::writeStats(Info, addFeMesh);
576 
577  feMesh.add(addFeMesh);
578  }
579 
580  if (selfIntersect != surfaceIntersection::NONE)
581  {
582  triSurfaceSearch query(surf);
583  surfaceIntersection intersect(query, surfaceDict);
584 
585  // Remove rounding noise. For consistency could use 1e-6,
586  // as per extractFromFile implementation
587 
588  intersect.mergePoints(10*SMALL);
589 
590  labelPair sizeInfo
591  (
592  intersect.cutPoints().size(),
593  intersect.cutEdges().size()
594  );
595 
596  if (intersect.cutEdges().size())
597  {
598  extendedEdgeMesh addMesh
599  (
600  intersect.cutPoints(),
601  intersect.cutEdges()
602  );
603 
604  feMesh.add(addMesh);
605 
606  sizeInfo[0] = addMesh.points().size();
607  sizeInfo[1] = addMesh.edges().size();
608  }
609  Info<< nl
610  << "intersection: "
612  << nl
613  << " points : " << sizeInfo[0] << nl
614  << " edges : " << sizeInfo[1] << nl;
615  }
616 
617  Info<< nl << "Final ";
619 
620  Info<< nl << "Writing extendedFeatureEdgeMesh to "
621  << feMesh.objectPath() << endl;
622 
623  mkDir(feMesh.path());
624 
625  if (writeObj)
626  {
627  feMesh.writeObj(feMesh.path()/outputName);
628  }
629 
630  feMesh.write();
631 
632  // Write a featureEdgeMesh (.eMesh) for backwards compatibility
633  // Used by snappyHexMesh (JUN-2017)
634  if (true)
635  {
636  featureEdgeMesh bfeMesh
637  (
638  IOobject
639  (
640  outputName + ".eMesh", // name
641  runTime.constant(), // instance
642  "triSurface",
643  runTime, // registry
646  false
647  ),
648  feMesh.points(),
649  feMesh.edges()
650  );
651 
652  Info<< nl << "Writing featureEdgeMesh to "
653  << bfeMesh.objectPath() << endl;
654 
655  bfeMesh.regIOobject::write();
656  }
657 
658 
659  // Output information
660 
661  const bool optCloseness =
662  surfaceDict.getOrDefault("closeness", false);
663 
664  const bool optProximity =
665  surfaceDict.getOrDefault("featureProximity", false);
666 
667  const bool optCurvature =
668  surfaceDict.getOrDefault("curvature", false);
669 
670 
671  // For VTK legacy format, we would need an a priori count of
672  // CellData and PointData fields.
673  // For convenience, we therefore only use the XML formats
674 
675  autoPtr<vtk::surfaceWriter> vtkWriter;
676 
677  if (optCloseness || optProximity || optCurvature)
678  {
679  if (writeVTK)
680  {
681  vtkWriter.reset
682  (
684  (
685  surf.points(),
686  faces,
687  (vtkOutputDir / outputName),
688  false // serial only
689  )
690  );
691 
692  vtkWriter->writeGeometry();
693 
694  Info<< "Writing VTK to "
695  << runTime.relativePath(vtkWriter->output()) << nl;
696  }
697  }
698  else
699  {
700  continue; // Nothing to output
701  }
702 
703 
704  // Option: "closeness"
705  if (optCloseness)
706  {
707  Pair<tmp<scalarField>> tcloseness =
709  (
710  runTime,
711  outputName,
712  surf,
713  45, // internalAngleTolerance
714  10 // externalAngleTolerance
715  );
716 
717  if (vtkWriter)
718  {
719  vtkWriter->beginCellData();
720  vtkWriter->write("internalCloseness", tcloseness[0]());
721  vtkWriter->write("externalCloseness", tcloseness[1]());
722  }
723  }
724 
725  // Option: "featureProximity"
726  if (optCloseness)
727  {
728  const scalar maxProximity =
729  surfaceDict.getOrDefault<scalar>("maxFeatureProximity", 1);
730 
731  tmp<scalarField> tproximity =
733  (
734  runTime,
735  outputName,
736  feMesh,
737  surf,
738  maxProximity
739  );
740 
741  if (vtkWriter)
742  {
743  vtkWriter->beginCellData();
744  vtkWriter->write("featureProximity", tproximity());
745  }
746  }
747 
748  // Option: "curvature"
749  if (optCurvature)
750  {
751  tmp<scalarField> tcurvature =
753  (
754  runTime,
755  outputName,
756  surf
757  );
758 
759  if (vtkWriter)
760  {
761  vtkWriter->beginPointData();
762  vtkWriter->write("curvature", tcurvature());
763  }
764  }
765 
766  Info<< endl;
767  }
768 
770 
771  Info<< "End\n" << endl;
772 
773  return 0;
774 }
775 
776 
777 // ************************************************************************* //
Foam::triSurface::writeStats
void writeStats(Ostream &os) const
Definition: triSurfaceIO.C:346
Foam::PrimitivePatch::points
const Field< point_type > & points() const noexcept
Definition: PrimitivePatch.H:295
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:63
Foam::word::lessExt
word lessExt() const
Definition: word.C:106
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
meshTools.H
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:190
Foam::triSurfaceLoader::loadingOption
loadingOption
Definition: triSurfaceLoader.H:63
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::fileName
A class for handling file names.
Definition: fileName.H:71
Foam::OBJstream
OFstream that keeps track of vertices.
Definition: OBJstream.H:54
setSystemRunTimeDictionaryIO.H
Foam::OFstream::name
virtual const fileName & name() const
Definition: OSstream.H:103
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryI.H:80
surfaceFeaturesExtraction.H
point.H
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:82
dictName
const word dictName("faMeshDefinition")
Foam::triSurfaceLoader::loadingOptionNames
static const Enum< loadingOption > loadingOptionNames
Definition: triSurfaceLoader.H:72
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
unitConversion.H
Unit conversion functions.
Foam::surfaceIntersection::intersectionType
intersectionType
Definition: surfaceIntersection.H:135
Foam::triSurfaceTools::writeCurvature
static tmp< scalarField > writeCurvature(const Time &runTime, const word &basename, const triSurface &surf)
Definition: triSurfaceCurvature.C:342
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:100
foamVtkSurfaceWriter.H
triSurface.H
Foam::OBJstream::write
virtual Ostream & write(const char c)
Definition: OBJstream.C:71
Foam::triSurface::patches
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:395
Foam::triSurfaceSearch
Helper class to search on triSurface.
Definition: triSurfaceSearch.H:54
Foam::TimePaths::constantPath
fileName constantPath() const
Definition: TimePathsI.H:123
Foam::autoPtr::empty
bool empty() const noexcept
Definition: autoPtr.H:274
Foam::surfaceFeaturesExtraction::method::New
static autoPtr< method > New(const dictionary &dict)
extendedFeatureEdgeMesh.H
Foam::triSurfaceLoader::FILE_REGION
@ FILE_REGION
"file" = One region for each file
Definition: triSurfaceLoader.H:66
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::writeVTK
void writeVTK(OFstream &os, const Type &value)
surfaceIntersection.H
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:85
Foam::edgeMeshTools::writeFeatureProximity
tmp< scalarField > writeFeatureProximity(const Time &runTime, const word &basename, const extendedEdgeMesh &emesh, const triSurface &surf, const scalar searchDistance)
Definition: edgeMeshFeatureProximity.C:184
outputName
word outputName("finiteArea-edges.obj")
Foam::Time::printExecutionTime
Ostream & printExecutionTime(OSstream &os) const
Definition: TimeIO.C:611
Foam::TimePaths::relativePath
fileName relativePath(const fileName &input, const bool caseTag=false) const
Definition: TimePathsI.H:80
Foam::argList::noFunctionObjects
static void noFunctionObjects(bool addWithOption=false)
Definition: argList.C:466
Foam::word::hasExt
bool hasExt() const
Definition: stringI.H:49
Foam::surfaceIntersection
Basic surface-surface intersection description. Constructed from two surfaces it creates a descriptio...
Definition: surfaceIntersection.H:130
plane.H
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:72
treeBoundBox.H
Foam::Info
messageStream Info
Foam::List::setSize
void setSize(const label n)
Definition: List.H:218
argList.H
Foam::Enum::getOrDefault
EnumType getOrDefault(const word &key, const dictionary &dict, const EnumType deflt, const bool failsafe=false) const
Definition: Enum.C:172
Foam::surfaceIntersection::NONE
@ NONE
None = invalid (for input only)
Definition: surfaceIntersection.H:141
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:453
Foam::surfaceFeatures
Holds feature edges/points of surface.
Definition: surfaceFeatures.H:69
patchNames
wordList patchNames(nPatches)
Foam::edgeMeshTools::writeStats
void writeStats(Ostream &os, const extendedFeatureEdgeMesh &emesh)
Definition: edgeMeshTools.C:29
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:379
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
os
OBJstream os(runTime.globalPath()/outputName)
dictIO
IOobject dictIO
Definition: setConstantMeshDictionaryIO.H:1
Foam
Definition: atmBoundaryLayer.C:26
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Definition: FlatOutput.H:217
Foam::triSurfaceLoader
Convenience class for loading single or multiple surface files from the constant/triSurface (or other...
Definition: triSurfaceLoader.H:58
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::dictionary::isDict
bool isDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryI.H:140
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Foam::IOobject::name
const word & name() const noexcept
Definition: IOobjectI.H:58
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
setRootCase.H
Foam::extendedEdgeMesh
Description of feature edges and points.
Definition: extendedEdgeMesh.H:81
Foam::surfaceIntersection::SELF_REGION
@ SELF_REGION
Self-intersection, region-wise only.
Definition: surfaceIntersection.H:140
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::triSurfaceTools::writeCloseness
static Pair< tmp< scalarField > > writeCloseness(const Time &runTime, const word &basename, const triSurface &surf, const scalar internalAngleTolerance=45, const scalar externalAngleTolerance=10)
Definition: triSurfaceCloseness.C:89
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:50
Foam::vtk::surfaceWriter
Write faces/points (optionally with fields) as a vtp file or a legacy vtk file.
Definition: foamVtkSurfaceWriter.H:61
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:47
Foam::Switch::OFF
@ OFF
Definition: Switch.H:92
createTime.H
edgeMeshTools.H
Foam::extendedFeatureEdgeMesh
extendedEdgeMesh + IO.
Definition: extendedFeatureEdgeMesh.H:52
patchId
label patchId(-1)
featureEdgeMesh.H
Foam::argList::noParallel
static void noParallel()
Definition: argList.C:503
Foam::triSurfaceLoader::OFFSET_REGION
@ OFFSET_REGION
"offset" = Offset regions per file
Definition: triSurfaceLoader.H:67
Foam::TimePaths::constant
const word & constant() const
Definition: TimePathsI.H:89
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:141
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
Foam::argList::addOption
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Definition: argList.C:328
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Definition: POSIX.C:533
Foam::featureEdgeMesh
edgeMesh + IO.
Definition: featureEdgeMesh.H:50
triSurfaceMesh.H
triSurfaceTools.H
triSurfaceLoader.H
Foam::objectRegistry::time
const Time & time() const noexcept
Definition: objectRegistry.H:174
OBJstream.H
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:398
Foam::surfaceIntersection::selfIntersectionNames
static const Enum< intersectionType > selfIntersectionNames
Definition: surfaceIntersection.H:145
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181