Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Application
28 surfaceFeatureExtract
29
30Group
31 grpSurfaceUtilities
32
33Description
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
181Note
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"
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
205using namespace Foam;
206
207// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208
209int 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
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.has_ext())
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 const dictionary* subDictPtr = nullptr;
378
379 //
380 // Extract features using the preferred extraction method
381 //
382 autoPtr<surfaceFeatures> features = extractor().features(surf);
383
384 // Trim set
385 // ~~~~~~~~
386
387 // Option: "trimFeatures" (dictionary)
388 if ((subDictPtr = surfaceDict.findDict("trimFeatures")) != nullptr)
389 {
390 const dictionary& trimDict = *subDictPtr;
391
392 const scalar minLen =
393 trimDict.getOrDefault<scalar>("minLen", 0);
394 const label minElem =
395 trimDict.getOrDefault<label>("minElem", 0);
396
397 // Trim away small groups of features
398 if (minLen > 0 || minElem > 0)
399 {
400 if (minLen > 0)
401 {
402 Info<< "Removing features of length < "
403 << minLen << endl;
404 }
405 if (minElem > 0)
406 {
407 Info<< "Removing features with number of edges < "
408 << minElem << endl;
409 }
410
411 features().trimFeatures
412 (
413 minLen, minElem, extractor().includedAngle()
414 );
415 }
416 }
417
418 // Subset
419 // ~~~~~~
420
421 // Convert to marked edges, points
422 List<surfaceFeatures::edgeStatus> edgeStat(features().toStatus());
423
424 // Option: "subsetFeatures" (dictionary)
425 if ((subDictPtr = surfaceDict.findDict("subsetFeatures")) != nullptr)
426 {
427 const dictionary& subsetDict = *subDictPtr;
428
429 treeBoundBox bb;
430
431 // Suboption: "insideBox"
432 if (subsetDict.readIfPresent("insideBox", bb))
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.readIfPresent("outsideBox", bb))
449 {
450 Info<< "Exclude edges outside box " << bb << endl;
451 features().excludeBox(edgeStat, bb);
452
453 {
454 OBJstream os("deleteBox.obj");
455
456 Info<< "Dumping bounding box " << bb
457 << " as lines to obj file "
458 << os.name() << endl;
459
460 os.write(bb);
461 }
462 }
463
464 // Suboption: "nonManifoldEdges" (false: remove non-manifold edges)
465 if (!subsetDict.getOrDefault("nonManifoldEdges", true))
466 {
467 Info<< "Removing all non-manifold edges"
468 << " (edges with > 2 connected faces) unless they"
469 << " cross multiple regions" << endl;
470
471 features().checkFlatRegionEdge
472 (
473 edgeStat,
474 1e-5, // tol
475 extractor().includedAngle()
476 );
477 }
478
479 if (!subsetDict.getOrDefault("strictNonManifoldEdges", true))
480 {
481 Info<< "Removing all non-manifold edges"
482 << " (edges with > 2 connected faces)"
483 << endl;
484
485 features().excludeNonManifold(edgeStat);
486 }
487
488 // Suboption: "openEdges" (false: remove open edges)
489 if (!subsetDict.getOrDefault("openEdges", true))
490 {
491 Info<< "Removing all open edges"
492 << " (edges with 1 connected face)" << endl;
493
494 features().excludeOpen(edgeStat);
495 }
496
497 // Suboption: "plane"
498 if (subsetDict.found("plane"))
499 {
500 plane cutPlane(subsetDict.lookup("plane"));
501
502 Info<< "Only include feature edges that intersect the plane"
503 << " with normal " << cutPlane.normal()
504 << " and origin " << cutPlane.origin() << endl;
505
506 features().subsetPlane(edgeStat, cutPlane);
507 }
508 }
509
510 surfaceFeatures newSet(surf);
511 newSet.setFromStatus(edgeStat, extractor().includedAngle());
512
513 Info<< nl << "Initial ";
514 newSet.writeStats(Info);
515
516 boolList surfBaffleRegions(surf.patches().size(), false);
517 if (surfaceDict.found("baffles"))
518 {
519 wordRes baffleSelect(surfaceDict.get<wordRes>("baffles"));
520
522 forAll(surf.patches(), patchi)
523 {
524 patchNames[patchi] = surf.patches()[patchi].name();
525 }
526
527 labelList indices(baffleSelect.matching(patchNames));
528
529 for (const label patchId : indices)
530 {
531 surfBaffleRegions[patchId] = true;
532 }
533
534 if (indices.size())
535 {
536 Info<< "Adding " << indices.size() << " baffle regions: (";
537
538 forAll(surfBaffleRegions, patchi)
539 {
540 if (surfBaffleRegions[patchi])
541 {
542 Info<< ' ' << patchNames[patchi];
543 }
544 }
545 Info<< " )" << nl << nl;
546 }
547 }
548
549 // Extracting and writing a extendedFeatureEdgeMesh
551 (
552 newSet,
553 runTime,
554 outputName + ".extendedFeatureEdgeMesh",
555 surfBaffleRegions
556 );
557
558
559 if ((subDictPtr = surfaceDict.findDict("addFeatures")) != nullptr)
560 {
561 const dictionary& addFeaturesDict = *subDictPtr;
562
563 const word addFeName = addFeaturesDict.get<word>("name");
564
565 Info<< "Adding (without merging) features from " << addFeName
566 << nl << endl;
567
569 (
571 (
572 addFeName,
573 runTime.time().constant(),
574 "extendedFeatureEdgeMesh",
575 runTime.time(),
578 )
579 );
580 Info<< "Read " << addFeMesh.name() << nl;
582
583 feMesh.add(addFeMesh);
584 }
585
586 if (selfIntersect != surfaceIntersection::NONE)
587 {
588 triSurfaceSearch query(surf);
589 surfaceIntersection intersect(query, surfaceDict);
590
591 // Remove rounding noise. For consistency could use 1e-6,
592 // as per extractFromFile implementation
593
594 intersect.mergePoints(10*SMALL);
595
596 labelPair sizeInfo
597 (
598 intersect.cutPoints().size(),
599 intersect.cutEdges().size()
600 );
601
602 if (intersect.cutEdges().size())
603 {
604 extendedEdgeMesh addMesh
605 (
606 intersect.cutPoints(),
607 intersect.cutEdges()
608 );
609
610 feMesh.add(addMesh);
611
612 sizeInfo[0] = addMesh.points().size();
613 sizeInfo[1] = addMesh.edges().size();
614 }
615 Info<< nl
616 << "intersection: "
618 << nl
619 << " points : " << sizeInfo[0] << nl
620 << " edges : " << sizeInfo[1] << nl;
621 }
622
623 Info<< nl << "Final ";
625
626 Info<< nl << "Writing extendedFeatureEdgeMesh to "
627 << feMesh.objectPath() << endl;
628
629 mkDir(feMesh.path());
630
631 if (writeObj)
632 {
633 feMesh.writeObj(feMesh.path()/outputName);
634 }
635
636 feMesh.write();
637
638 // Write a featureEdgeMesh (.eMesh) for backwards compatibility
639 // Used by snappyHexMesh (JUN-2017)
640 if (true)
641 {
642 featureEdgeMesh bfeMesh
643 (
645 (
646 outputName + ".eMesh", // name
647 runTime.constant(), // instance
648 "triSurface",
649 runTime, // registry
653 ),
654 feMesh.points(),
655 feMesh.edges()
656 );
657
658 Info<< nl << "Writing featureEdgeMesh to "
659 << bfeMesh.objectPath() << endl;
660
661 bfeMesh.regIOobject::write();
662 }
663
664
665 // Output information
666
667 const bool optCloseness =
668 surfaceDict.getOrDefault("closeness", false);
669
670 const bool optProximity =
671 surfaceDict.getOrDefault("featureProximity", false);
672
673 const bool optCurvature =
674 surfaceDict.getOrDefault("curvature", false);
675
676
677 // For VTK legacy format, we would need an a priori count of
678 // CellData and PointData fields.
679 // For convenience, we therefore only use the XML formats
680
682
683 if (optCloseness || optProximity || optCurvature)
684 {
685 if (writeVTK)
686 {
687 vtkWriter.reset
688 (
690 (
691 surf.points(),
692 faces,
693 (vtkOutputDir / outputName),
694 false // serial only
695 )
696 );
697
698 vtkWriter->writeGeometry();
699
700 Info<< "Writing VTK to "
701 << runTime.relativePath(vtkWriter->output()) << nl;
702 }
703 }
704 else
705 {
706 continue; // Nothing to output
707 }
708
709
710 // Option: "closeness"
711 if (optCloseness)
712 {
713 Pair<tmp<scalarField>> tcloseness =
715 (
716 runTime,
718 surf,
719 45, // internalAngleTolerance
720 10 // externalAngleTolerance
721 );
722
723 if (vtkWriter)
724 {
725 vtkWriter->beginCellData();
726 vtkWriter->write("internalCloseness", tcloseness[0]());
727 vtkWriter->write("externalCloseness", tcloseness[1]());
728 }
729 }
730
731 // Option: "featureProximity"
732 if (optCloseness)
733 {
734 const scalar maxProximity =
735 surfaceDict.getOrDefault<scalar>("maxFeatureProximity", 1);
736
737 tmp<scalarField> tproximity =
739 (
740 runTime,
742 feMesh,
743 surf,
744 maxProximity
745 );
746
747 if (vtkWriter)
748 {
749 vtkWriter->beginCellData();
750 vtkWriter->write("featureProximity", tproximity());
751 }
752 }
753
754 // Option: "curvature"
755 if (optCurvature)
756 {
757 tmp<scalarField> tcurvature =
759 (
760 runTime,
762 surf
763 );
764
765 if (vtkWriter)
766 {
767 vtkWriter->beginPointData();
768 vtkWriter->write("curvature", tcurvature());
769 }
770 }
771
772 Info<< endl;
773 }
774
775 runTime.printExecutionTime(Info);
776
777 Info<< "End\n" << endl;
778
779 return 0;
780}
781
782
783// ************************************************************************* //
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
void setSize(label n)
Alias for resize().
Definition List.H:536
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
const Field< point_type > & points() const noexcept
Return reference to global points.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
static void noParallel()
Remove the parallel options.
Definition argList.C:599
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream.
Definition dictionary.C:367
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Description of feature edges and points.
A class for handling file names.
Definition fileName.H:75
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
static autoPtr< method > New(const dictionary &dict)
Select constructed from dictionary.
Holds feature edges/points of surface.
Basic surface-surface intersection description. Constructed from two surfaces it creates a descriptio...
intersectionType
Surface intersection types for classify, doCutEdges.
@ SELF_REGION
Self-intersection, region-wise only.
@ NONE
None = invalid (for input only).
static const Enum< intersectionType > selfIntersectionNames
The user-selectable self-intersection enumeration names.
A class for managing temporary objects.
Definition tmp.H:75
Standard boundBox with extra functionality for use in octree.
Convenience class for loading single or multiple surface files from the constant/triSurface (or other...
static const Enum< loadingOption > loadingOptionNames
The loading enumeration names.
loadingOption
The file loading options for triSurfaceLoader.
@ FILE_REGION
"file" = One region for each file
@ OFFSET_REGION
"offset" = Offset regions per file
Helper class to search on triSurface.
static tmp< scalarField > writeCurvature(const Time &runTime, const word &basename, const triSurface &surf)
Write surface curvature at the vertex points and return the field.
static Pair< tmp< scalarField > > writeCloseness(const Time &runTime, const word &basename, const triSurface &surf, const scalar internalAngleTolerance=45, const scalar externalAngleTolerance=10)
Check and write internal/external closeness fields.
Triangulated surface description with patch information.
Definition triSurface.H:74
const geometricSurfacePatchList & patches() const noexcept
Definition triSurface.H:509
void writeStats(Ostream &os) const
Write some statistics.
Write faces/points (optionally with fields) as a vtp file or a legacy vtk file.
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
word outputName("finiteArea-edges.obj")
OBJstream os(runTime.globalPath()/outputName)
const word dictName("faMeshDefinition")
label patchId(-1)
void writeStats(Ostream &os, const extendedFeatureEdgeMesh &emesh)
Write some information.
tmp< scalarField > writeFeatureProximity(const Time &runTime, const word &basename, const extendedEdgeMesh &emesh, const triSurface &surf, const scalar searchDistance)
Calculate proximity of the features to the surface and write the field.
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
void writeVTK(OFstream &os, const Type &value)
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
List< bool > boolList
A List of bools.
Definition List.H:60
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
dictionary dict
IOobject dictIO
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.
mkDir(pdfPath)