STLReader.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) 2016-2020 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
27\*---------------------------------------------------------------------------*/
28
29#include "STLReader.H"
30#include "STLAsciiParse.H"
31#include "Map.H"
32#include "IFstream.H"
33#include "mergePoints.H"
34
35#undef DEBUG_STLBINARY
36
37/* * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * */
38
40(
41 Foam::debug::optimisationSwitch("fileFormats::stl", 0)
42);
43
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
48(
49 Detail::STLAsciiParse& parsed
50)
51{
52 sorted_ = parsed.sorted();
53
54 points_.transfer(parsed.points());
55 zoneIds_.transfer(parsed.facets());
56 names_.transfer(parsed.names());
57 sizes_.transfer(parsed.sizes());
58
59 format_ = STLFormat::ASCII;
60
61 parsed.clear();
62}
63
64
65bool Foam::fileFormats::STLReader::readASCII
66(
67 const fileName& filename
68)
69{
70 // No runtime selection of parser (only via optimisationSwitch)
71 // this is something that is infrequently changed.
72 if (parserType == 1)
73 {
74 return readAsciiRagel(filename);
75 }
76 else if (parserType == 2)
77 {
78 return readAsciiManual(filename);
79 }
80 return readAsciiFlex(filename);
81}
82
83
84bool Foam::fileFormats::STLReader::readBINARY
85(
86 const fileName& filename
87)
88{
89 sorted_ = true;
90 format_ = STLFormat::UNKNOWN;
91
92 label nTris = 0;
93 std::unique_ptr<std::istream> streamPtr
94 {
95 readBinaryHeader(filename, nTris)
96 };
97
98 if (!streamPtr)
99 {
101 << "Error reading file " << filename
102 << " or file " << filename + ".gz"
103 << exit(FatalError);
104 }
105 auto& is = *streamPtr;
106
107#ifdef DEBUG_STLBINARY
108 Info<< "# " << nTris << " facets" << endl;
109 label prevZone = -1;
110#endif
111
112 points_.setSize(3*nTris);
113 zoneIds_.setSize(nTris);
114
115 Map<label> lookup;
116 DynamicList<label> dynSizes;
117
118 label ptI = 0;
119 label zoneI = -1;
120 forAll(zoneIds_, facei)
121 {
122 // Read STL triangle
123 STLtriangle stlTri(is);
124
125 // transcribe the vertices of the STL triangle -> points
126 points_[ptI++] = stlTri.a();
127 points_[ptI++] = stlTri.b();
128 points_[ptI++] = stlTri.c();
129
130 // interpret STL attribute as a zone
131 const label origId = stlTri.attrib();
132
133 auto fnd = lookup.cfind(origId);
134 if (fnd.found())
135 {
136 if (zoneI != *fnd)
137 {
138 sorted_ = false; // Group appeared out of order
139 }
140 zoneI = *fnd;
141 }
142 else
143 {
144 zoneI = dynSizes.size();
145 lookup.insert(origId, zoneI);
146 dynSizes.append(0);
147 }
148
149 zoneIds_[facei] = zoneI;
150 dynSizes[zoneI]++;
151
152#ifdef DEBUG_STLBINARY
153 if (prevZone != zoneI)
154 {
155 if (prevZone != -1)
156 {
157 Info<< "endsolid zone" << prevZone << nl;
158 }
159 prevZone = zoneI;
160
161 Info<< "solid zone" << prevZone << nl;
162 }
163
164 stlTri.print(Info);
165#endif
166 }
167
168#ifdef DEBUG_STLBINARY
169 if (prevZone != -1)
170 {
171 Info<< "endsolid zone" << prevZone << nl;
172 }
173#endif
174
175 names_.clear();
176 sizes_.transfer(dynSizes);
177
178 format_ = STLFormat::BINARY;
179 return true;
180}
181
182
183bool Foam::fileFormats::STLReader::readFile
184(
185 const fileName& filename,
186 const STLFormat format
187)
188{
189 if
190 (
191 format == STLFormat::UNKNOWN
192 ? detectBinaryHeader(filename)
193 : format == STLFormat::BINARY
194 )
195 {
196 return readBINARY(filename);
197 }
198 else
199 {
200 return readASCII(filename);
201 }
202}
203
204
205// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
206
208(
209 const fileName& filename
210)
211:
212 sorted_(true),
213 points_(),
214 zoneIds_(),
215 names_(),
216 sizes_(),
217 format_(STLFormat::UNKNOWN)
218{
219 // Auto-detect ASCII/BINARY format
220 readFile(filename, STLFormat::UNKNOWN);
221}
222
223
225(
226 const fileName& filename,
227 const STLFormat format
228)
229:
230 sorted_(true),
231 points_(),
232 zoneIds_(),
233 names_(),
234 sizes_(),
235 format_(STLFormat::UNKNOWN)
236{
237 // Manually specified ASCII/BINARY format
238 readFile(filename, format);
239}
240
241
242// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
243
245{
246 sorted_ = true;
247 points_.clear();
248 zoneIds_.clear();
249 names_.clear();
250 sizes_.clear();
251 format_ = STLFormat::UNKNOWN;
252}
253
254
256(
257 labelList& pointMap
258) const
259{
260 // With the merge distance depending on the input format (ASCII | BINARY),
261 // but must be independent of WM_SP or WM_DP flag.
262 // - floatScalarSMALL = 1e-6
263 // - doubleScalarSMALL = 1e-15
264
265 return mergePointsMap
266 (
267 (format_ == STLFormat::BINARY ? 10 : 100) * doubleScalarSMALL,
268 pointMap
269 );
270}
271
272
274(
275 const scalar mergeTol,
276 labelList& pointMap
277) const
278{
279 return Foam::mergePoints
280 (
281 points_,
282 mergeTol,
283 false, // verbose
284 pointMap
285 );
286}
287
288
289// ************************************************************************* //
void transfer(List< T > &list)
Definition: List.C:447
STLFormat
Enumeration for the format of data in the stream.
Definition: STLCore.H:62
@ UNKNOWN
Detect based on (input) content or (output) extension.
Definition: STLCore.H:65
Internal class used by the STLsurfaceFormat and triSurface.
Definition: STLReader.H:70
label mergePointsMap(labelList &pointMap) const
Calculate merge points mapping, return old to new pointMap.
Definition: STLReader.C:256
static int parserType
ASCII parser types (0=Flex, 1=Ragel, 2=Manual)
Definition: STLReader.H:127
void clear()
Flush all values.
Definition: STLReader.C:244
A class for handling file names.
Definition: fileName.H:76
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Geometric merging of points. See below.
int optimisationSwitch(const char *name, const int deflt=0)
Lookup optimisation switch or add default value.
Definition: debug.C:237
messageStream Info
Information stream (stdout output on master, null elsewhere)
label mergePoints(const PointList &points, labelList &pointToUnique, labelList &uniquePoints, const scalar mergeTol=SMALL, const bool verbose=false)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
error FatalError
constexpr doubleScalar doubleScalarSMALL
Definition: doubleScalar.H:61
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
word format(conversionProperties.get< word >("format"))
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333