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