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  // Clear everything
57  this->clear();
58  sorted_ = true;
59 
60  IFstream is(filename);
61  if (!is.good())
62  {
64  << "Cannot read file " << filename << nl
65  << exit(FatalError);
66  }
67 
68  // uses similar structure as STL, just some points
69  // the rest of the reader resembles the STL binary reader
70  DynamicList<STLpoint> dynPoints;
71  DynamicList<label> dynZones;
72  DynamicList<word> dynNames;
73  DynamicList<label> dynSizes;
74  HashTable<label> lookup;
75 
76  // place faces without a group in zone0
77  label zoneI = 0;
78  dynSizes.append(zoneI);
79  lookup.insert("zoneI", zoneI);
80 
81  while (is.good())
82  {
83  string line = this->getLineNoComment(is);
84 
85  if (line.empty())
86  {
87  break;
88  }
89 
90  // Do not handle continuations
91 
92  IStringStream lineStream(line);
93 
94  STLpoint p(getSTLpoint(lineStream));
95 
96  if (!lineStream) break;
97 
98  dynPoints.append(p);
99  dynPoints.append(getSTLpoint(lineStream));
100  dynPoints.append(getSTLpoint(lineStream));
101 
102  // zone/colour in .tri file starts with 0x. Skip.
103  // ie, instead of having 0xFF, skip 0 and leave xFF to
104  // get read as a word and name it "zoneFF"
105 
106  char zeroChar;
107  lineStream >> zeroChar;
108 
109  const word rawName(lineStream);
110  const word name("zone" + rawName.substr(1));
111 
112  const auto iter = lookup.cfind(name);
113  if (iter.found())
114  {
115  if (zoneI != iter.val())
116  {
117  sorted_ = false; // Group appeared out of order
118  zoneI = iter.val();
119  }
120  }
121  else
122  {
123  zoneI = dynSizes.size();
124  if (lookup.insert(name, zoneI))
125  {
126  dynNames.append(name);
127  dynSizes.append(0);
128  }
129  }
130 
131  dynZones.append(zoneI);
132  dynSizes[zoneI]++;
133  }
134 
135  // skip empty groups
136  label nZone = 0;
137  forAll(dynSizes, zonei)
138  {
139  if (dynSizes[zonei])
140  {
141  if (nZone != zonei)
142  {
143  dynNames[nZone] = dynNames[zonei];
144  dynSizes[nZone] = dynSizes[zonei];
145  }
146  ++nZone;
147  }
148  }
149 
150  // truncate addressed size
151  dynNames.setCapacity(nZone);
152  dynSizes.setCapacity(nZone);
153 
154  // transfer to normal lists
155  points_.transfer(dynPoints);
156  zoneIds_.transfer(dynZones);
157  names_.transfer(dynNames);
158  sizes_.transfer(dynSizes);
159 
160  return true;
161 }
162 
163 
164 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
165 
167 (
168  const fileName& filename
169 )
170 :
171  sorted_(true),
172  points_(),
173  zoneIds_(),
174  names_(),
175  sizes_()
176 {
177  readFile(filename);
178 }
179 
180 
181 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
182 
184 {
185  sorted_ = true;
186  points_.clear();
187  zoneIds_.clear();
188  names_.clear();
189  sizes_.clear();
190 }
191 
192 
194 (
195  labelList& pointMap
196 ) const
197 {
198  // Use merge tolerance as per STL ASCII
199  return mergePointsMap
200  (
201  100 * doubleScalarSMALL,
202  pointMap
203  );
204 }
205 
206 
208 (
209  const scalar mergeTol,
210  labelList& pointMap
211 ) const
212 {
213  return Foam::mergePoints
214  (
215  points_,
216  mergeTol,
217  false, // verbose
218  pointMap
219  );
220 }
221 
222 
223 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
TRIReader.H
Foam::doubleScalarSMALL
constexpr doubleScalar doubleScalarSMALL
Definition: doubleScalar.H:61
StringStream.H
Input/output from string buffers.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Map.H
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::fileFormats::surfaceFormatsCore::getLineNoComment
static string getLineNoComment(ISstream &is, const char comment='#')
Read non-empty and non-comment line.
Definition: surfaceFormatsCore.C:43
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
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:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::fileFormats::TRIReader::clear
void clear()
Flush all values.
Definition: TRIReader.C:183
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:194
Foam::fileFormats::TRIReader::TRIReader
TRIReader(const fileName &filename)
Read from file, filling in the information.
Definition: TRIReader.C:167
mergePoints.H
Merge points. See below.
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::getSTLpoint
static STLpoint getSTLpoint(Istream &is)
Definition: TRIReader.C:41
surfaceFormatsCore.H