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-2017 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 \*---------------------------------------------------------------------------*/
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 
47 void Foam::fileFormats::STLReader::transfer
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 
65 bool 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 
84 bool Foam::fileFormats::STLReader::readBINARY
85 (
86  const fileName& filename
87 )
88 {
89  sorted_ = true;
90  format_ = STLFormat::UNKNOWN;
91 
92  label nTris = 0;
93  autoPtr<istream> streamPtr = readBinaryHeader(filename, nTris);
94 
95  if (!streamPtr.valid())
96  {
98  << "Error reading file " << filename
99  << " or file " << filename + ".gz"
100  << exit(FatalError);
101  }
102 
103  istream& is = streamPtr();
104 
105 #ifdef DEBUG_STLBINARY
106  Info<< "# " << nTris << " facets" << endl;
107  label prevZone = -1;
108 #endif
109 
110  points_.setSize(3*nTris);
111  zoneIds_.setSize(nTris);
112 
113  Map<label> lookup;
114  DynamicList<label> dynSizes;
115 
116  label ptI = 0;
117  label zoneI = -1;
118  forAll(zoneIds_, facei)
119  {
120  // Read STL triangle
121  STLtriangle stlTri(is);
122 
123  // transcribe the vertices of the STL triangle -> points
124  points_[ptI++] = stlTri.a();
125  points_[ptI++] = stlTri.b();
126  points_[ptI++] = stlTri.c();
127 
128  // interpret STL attribute as a zone
129  const label origId = stlTri.attrib();
130 
131  Map<label>::const_iterator fnd = lookup.find(origId);
132  if (fnd != lookup.end())
133  {
134  if (zoneI != fnd())
135  {
136  // group appeared out of order
137  sorted_ = false;
138  }
139  zoneI = fnd();
140  }
141  else
142  {
143  zoneI = dynSizes.size();
144  lookup.insert(origId, zoneI);
145  dynSizes.append(0);
146  }
147 
148  zoneIds_[facei] = zoneI;
149  dynSizes[zoneI]++;
150 
151 #ifdef DEBUG_STLBINARY
152  if (prevZone != zoneI)
153  {
154  if (prevZone != -1)
155  {
156  Info<< "endsolid zone" << prevZone << nl;
157  }
158  prevZone = zoneI;
159 
160  Info<< "solid zone" << prevZone << nl;
161  }
162 
163  stlTri.print(Info);
164 #endif
165  }
166 
167 #ifdef DEBUG_STLBINARY
168  if (prevZone != -1)
169  {
170  Info<< "endsolid zone" << prevZone << nl;
171  }
172 #endif
173 
174  names_.clear();
175  sizes_.transfer(dynSizes);
176 
177  format_ = STLFormat::BINARY;
178  return true;
179 }
180 
181 
182 bool Foam::fileFormats::STLReader::readFile
183 (
184  const fileName& filename,
185  const STLFormat format
186 )
187 {
188  if
189  (
190  format == STLFormat::UNKNOWN
191  ? detectBinaryHeader(filename)
192  : format == STLFormat::BINARY
193  )
194  {
195  return readBINARY(filename);
196  }
197  else
198  {
199  return readASCII(filename);
200  }
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
205 
206 Foam::fileFormats::STLReader::STLReader
207 (
208  const fileName& filename
209 )
210 :
211  sorted_(true),
212  points_(),
213  zoneIds_(),
214  names_(),
215  sizes_(),
216  format_(STLFormat::UNKNOWN)
217 {
218  // Auto-detect ASCII/BINARY format
219  readFile(filename, STLFormat::UNKNOWN);
220 }
221 
222 
223 Foam::fileFormats::STLReader::STLReader
224 (
225  const fileName& filename,
226  const STLFormat format
227 )
228 :
229  sorted_(true),
230  points_(),
231  zoneIds_(),
232  names_(),
233  sizes_(),
234  format_(STLFormat::UNKNOWN)
235 {
236  // Manually specified ASCII/BINARY format
237  readFile(filename, format);
238 }
239 
240 
241 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
242 
244 {
245  sorted_ = true;
246  points_.clear();
247  zoneIds_.clear();
248  names_.clear();
249  sizes_.clear();
250  format_ = STLFormat::UNKNOWN;
251 }
252 
253 
255 (
256  labelList& pointMap
257 ) const
258 {
259  // With the merge distance depending on the input format (ASCII | BINARY),
260  // but must be independent of WM_SP or WM_DP flag.
261  // - floatScalarSMALL = 1e-6
262  // - doubleScalarSMALL = 1e-15
263 
264  return mergePointsMap
265  (
266  (format_ == STLFormat::BINARY ? 10 : 100) * doubleScalarSMALL,
267  pointMap
268  );
269 }
270 
271 
273 (
274  const scalar mergeTol,
275  labelList& pointMap
276 ) const
277 {
278  return Foam::mergePoints
279  (
280  points_,
281  mergeTol,
282  false, // verbose
283  pointMap
284  );
285 }
286 
287 
288 // ************************************************************************* //
Foam::Map< label >::const_iterator
typename parent_type::const_iterator const_iterator
Definition: Map.H:70
STLReader.H
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::doubleScalarSMALL
constexpr doubleScalar doubleScalarSMALL
Definition: doubleScalar.H:59
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Map.H
format
word format(conversionProperties.get< word >("format"))
Foam::label
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:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::cellModeller::lookup
const cellModel * lookup(const word &modelName)
Deprecated(2017-11) equivalent to cellModel::ptr static method.
Definition: cellModeller.H:46
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:436
Foam::fileFormats::STLReader::clear
void clear()
Flush all values.
Definition: STLReader.C:243
Foam::fileFormats::STLReader::mergePointsMap
label mergePointsMap(labelList &pointMap) const
Calculate merge points mapping, return old to new pointMap.
Definition: STLReader.C:255
IFstream.H
Foam::debug::optimisationSwitch
int optimisationSwitch(const char *name, const int deflt=0)
Lookup optimisation switch or add default value.
Definition: debug.C:237
Foam::FatalError
error FatalError
Foam::mergePoints
label mergePoints(const PointList &points, const scalar mergeTol, const bool verbose, labelList &pointMap, typename PointList::const_reference origin=PointList::value_type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::fileFormats::STLCore::STLFormat
STLFormat
Enumeration for the format of data in the stream.
Definition: STLCore.H:61
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::List< label >
mergePoints.H
Merge points. See below.
Foam::fileFormats::STLReader::parserType
static int parserType
ASCII parser types (0=Flex, 1=Ragel, 2=Manual)
Definition: STLReader.H:127
STLAsciiParse.H