TRIsurfaceFormatCore.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "TRIsurfaceFormat.H"
27 #include "IFstream.H"
28 #include "IOmanip.H"
29 #include "IStringStream.H"
30 #include "Map.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 Foam::fileFormats::TRIsurfaceFormatCore::TRIsurfaceFormatCore
35 (
36  const fileName& filename
37 )
38 :
39  sorted_(true),
40  points_(0),
41  zoneIds_(0),
42  sizes_(0)
43 {
44  read(filename);
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
49 
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 bool Foam::fileFormats::TRIsurfaceFormatCore::read
57 (
58  const fileName& filename
59 )
60 {
61  this->clear();
62  sorted_ = true;
63 
64  IFstream is(filename);
65  if (!is.good())
66  {
68  (
69  "fileFormats::TRIsurfaceFormatCore::read(const fileName&)"
70  )
71  << "Cannot read file " << filename
72  << exit(FatalError);
73  }
74 
75  // uses similar structure as STL, just some points
76  // the rest of the reader resembles the STL binary reader
77  DynamicList<point> dynPoints;
78  DynamicList<label> dynZones;
79  DynamicList<label> dynSizes;
81 
82  // place faces without a group in zone0
83  label zoneI = 0;
84  dynSizes.append(zoneI);
85  lookup.insert("zoneI", zoneI);
86 
87  while (is.good())
88  {
89  string line = this->getLineNoComment(is);
90 
91  // handle continuations ?
92  // if (line[line.size()-1] == '\\')
93  // {
94  // line.substr(0, line.size()-1);
95  // line += this->getLineNoComment(is);
96  // }
97 
98  IStringStream lineStream(line);
99 
100  point p
101  (
102  readScalar(lineStream),
103  readScalar(lineStream),
104  readScalar(lineStream)
105  );
106 
107  if (!lineStream) break;
108 
109  dynPoints.append(p);
110  dynPoints.append
111  (
112  point
113  (
114  readScalar(lineStream),
115  readScalar(lineStream),
116  readScalar(lineStream)
117  )
118  );
119  dynPoints.append
120  (
121  point
122  (
123  readScalar(lineStream),
124  readScalar(lineStream),
125  readScalar(lineStream)
126  )
127  );
128 
129  // zone/colour in .tri file starts with 0x. Skip.
130  // ie, instead of having 0xFF, skip 0 and leave xFF to
131  // get read as a word and name it "zoneFF"
132 
133  char zero;
134  lineStream >> zero;
135 
136  word rawName(lineStream);
137  word name("zone" + rawName(1, rawName.size()-1));
138 
139  HashTable<label>::const_iterator fnd = lookup.find(name);
140  if (fnd != lookup.end())
141  {
142  if (zoneI != fnd())
143  {
144  // group appeared out of order
145  sorted_ = false;
146  }
147  zoneI = fnd();
148  }
149  else
150  {
151  zoneI = dynSizes.size();
152  lookup.insert(name, zoneI);
153  dynSizes.append(0);
154  }
155 
156  dynZones.append(zoneI);
157  dynSizes[zoneI]++;
158  }
159 
160  // skip empty groups
161  label nZone = 0;
162  forAll(dynSizes, zoneI)
163  {
164  if (dynSizes[zoneI])
165  {
166  if (nZone != zoneI)
167  {
168  dynSizes[nZone] = dynSizes[zoneI];
169  }
170  nZone++;
171  }
172  }
173  // truncate addressed size
174  dynSizes.setCapacity(nZone);
175 
176  // transfer to normal lists
177  points_.transfer(dynPoints);
178  zoneIds_.transfer(dynZones);
179  sizes_.transfer(dynSizes);
180 
181  return true;
182 }
183 
184 
185 // ************************************************************************* //
Input from memory buffer stream.
Definition: IStringStream.H:49
void setCapacity(const label)
Alter the size of the underlying storage.
Definition: DynamicListI.H:109
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:390
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
A class for handling words, derived from string.
Definition: word.H:59
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:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:310
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Input from file stream.
Definition: IFstream.H:81
A line primitive.
Definition: line.H:56
stressControl lookup("compactNormalStress") >> compactNormalStress
volScalarField & p
Definition: createFields.H:51
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47
Istream and Ostream manipulators taking arguments.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
error FatalError
A class for handling file names.
Definition: fileName.H:69
static string getLineNoComment(IFstream &)
Read non-comment line.
bool read(const char *, int32_t &)
Definition: int32IO.C:87
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:106
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:139