TRIReader.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-2015 OpenFOAM Foundation
9 Copyright (C) 2017-2019 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 "TRIReader.H"
30#include "surfaceFormatsCore.H"
31#include "IFstream.H"
32#include "IOmanip.H"
33#include "StringStream.H"
34#include "DynamicList.H"
35#include "mergePoints.H"
36#include "Map.H"
37
38// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
39
40namespace Foam
41{
42static inline STLpoint getSTLpoint(Istream& is)
43{
44 scalar a = readScalar(is);
45 scalar b = readScalar(is);
46 scalar c = readScalar(is);
47
48 return STLpoint(a, b, c);
49}
50} // End namespace Foam
51
52
53// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54
55bool Foam::fileFormats::TRIReader::readFile(const fileName& filename)
56{
57 // Clear everything
58 this->clear();
59 sorted_ = true;
60
61 IFstream is(filename);
62 if (!is.good())
63 {
65 << "Cannot read file " << filename << nl
66 << exit(FatalError);
67 }
68
69 // uses similar structure as STL, just some points
70 // the rest of the reader resembles the STL binary reader
71 DynamicList<STLpoint> dynPoints;
72 DynamicList<label> dynZones;
73 DynamicList<word> dynNames;
74 DynamicList<label> dynSizes;
75 HashTable<label> lookup;
76
77 // place faces without a group in zone0
78 label zoneI = 0;
79 dynSizes.append(zoneI);
80 lookup.insert("zoneI", zoneI);
81
82 while (is.good())
83 {
84 string line = this->getLineNoComment(is);
85
86 if (line.empty())
87 {
88 break;
89 }
90
91 // Do not handle continuations
92
93 IStringStream lineStream(line);
94
95 STLpoint p(getSTLpoint(lineStream));
96
97 if (!lineStream) break;
98
99 dynPoints.append(p);
100 dynPoints.append(getSTLpoint(lineStream));
101 dynPoints.append(getSTLpoint(lineStream));
102
103 // zone/colour in .tri file starts with 0x. Skip.
104 // ie, instead of having 0xFF, skip 0 and leave xFF to
105 // get read as a word and name it "zoneFF"
106
107 char zeroChar;
108 lineStream >> zeroChar;
109
110 const word rawName(lineStream);
111 const word name("zone" + rawName.substr(1));
112
113 const auto iter = lookup.cfind(name);
114 if (iter.found())
115 {
116 if (zoneI != iter.val())
117 {
118 sorted_ = false; // Group appeared out of order
119 zoneI = iter.val();
120 }
121 }
122 else
123 {
124 zoneI = dynSizes.size();
125 if (lookup.insert(name, zoneI))
126 {
127 dynNames.append(name);
128 dynSizes.append(0);
129 }
130 }
131
132 dynZones.append(zoneI);
133 dynSizes[zoneI]++;
134 }
135
136 // skip empty groups
137 label nZone = 0;
138 forAll(dynSizes, zonei)
139 {
140 if (dynSizes[zonei])
141 {
142 if (nZone != zonei)
143 {
144 dynNames[nZone] = dynNames[zonei];
145 dynSizes[nZone] = dynSizes[zonei];
146 }
147 ++nZone;
148 }
149 }
150
151 // truncate addressed size
152 dynNames.setCapacity(nZone);
153 dynSizes.setCapacity(nZone);
154
155 // transfer to normal lists
156 points_.transfer(dynPoints);
157 zoneIds_.transfer(dynZones);
158 names_.transfer(dynNames);
159 sizes_.transfer(dynSizes);
160
161 return true;
162}
163
164
165// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
166
168(
169 const fileName& filename
170)
171:
172 sorted_(true),
173 points_(),
174 zoneIds_(),
175 names_(),
176 sizes_()
177{
178 readFile(filename);
179}
180
181
182// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
183
185{
186 sorted_ = true;
187 points_.clear();
188 zoneIds_.clear();
189 names_.clear();
190 sizes_.clear();
191}
192
193
195(
196 labelList& pointMap
197) const
198{
199 // Use merge tolerance as per STL ASCII
200 return mergePointsMap
201 (
202 100 * doubleScalarSMALL,
203 pointMap
204 );
205}
206
207
209(
210 const scalar mergeTol,
211 labelList& pointMap
212) const
213{
214 return Foam::mergePoints
215 (
216 points_,
217 mergeTol,
218 false, // verbose
219 pointMap
220 );
221}
222
223
224// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Input/output from string buffers.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
void transfer(List< T > &list)
Definition: List.C:447
A vertex point or facet normal representation for STL files.
Definition: STLpoint.H:52
TRI (triangle) file reader.
Definition: TRIReader.H:64
label mergePointsMap(labelList &pointMap) const
Calculate merge points mapping, return old to new pointMap.
Definition: TRIReader.C:195
void clear()
Flush all values.
Definition: TRIReader.C:184
static string getLineNoComment(ISstream &is, const char comment='#')
Read non-empty and non-comment line.
A class for handling file names.
Definition: fileName.H:76
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Geometric merging of points. See below.
Namespace for OpenFOAM.
static STLpoint getSTLpoint(Istream &is)
Definition: TRIReader.C:42
label mergePoints(const PointList &points, labelList &pointToUnique, labelList &uniquePoints, const scalar mergeTol=SMALL, const bool verbose=false)
error FatalError
constexpr doubleScalar doubleScalarSMALL
Definition: doubleScalar.H:61
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
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
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333