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 
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 
57 (
58  const fileName& filename
59 )
60 {
61  this->clear();
62  sorted_ = true;
63 
64  IFstream is(filename);
65  if (!is.good())
66  {
68  << "Cannot read file " << filename
69  << exit(FatalError);
70  }
71 
72  // uses similar structure as STL, just some points
73  // the rest of the reader resembles the STL binary reader
74  DynamicList<point> dynPoints;
75  DynamicList<label> dynZones;
76  DynamicList<label> dynSizes;
78 
79  // place faces without a group in zone0
80  label zoneI = 0;
81  dynSizes.append(zoneI);
82  lookup.insert("zoneI", zoneI);
83 
84  while (is.good())
85  {
86  string line = this->getLineNoComment(is);
87 
88  // handle continuations ?
89  // if (line[line.size()-1] == '\\')
90  // {
91  // line.substr(0, line.size()-1);
92  // line += this->getLineNoComment(is);
93  // }
94 
95  IStringStream lineStream(line);
96 
97  point p
98  (
99  readScalar(lineStream),
100  readScalar(lineStream),
101  readScalar(lineStream)
102  );
103 
104  if (!lineStream) break;
105 
106  dynPoints.append(p);
107  dynPoints.append
108  (
109  point
110  (
111  readScalar(lineStream),
112  readScalar(lineStream),
113  readScalar(lineStream)
114  )
115  );
116  dynPoints.append
117  (
118  point
119  (
120  readScalar(lineStream),
121  readScalar(lineStream),
122  readScalar(lineStream)
123  )
124  );
125 
126  // zone/colour in .tri file starts with 0x. Skip.
127  // ie, instead of having 0xFF, skip 0 and leave xFF to
128  // get read as a word and name it "zoneFF"
129 
130  char zero;
131  lineStream >> zero;
132 
133  word rawName(lineStream);
134  word name("zone" + rawName(1, rawName.size()-1));
135 
137  if (fnd != lookup.end())
138  {
139  if (zoneI != fnd())
140  {
141  // group appeared out of order
142  sorted_ = false;
143  }
144  zoneI = fnd();
145  }
146  else
147  {
148  zoneI = dynSizes.size();
149  lookup.insert(name, zoneI);
150  dynSizes.append(0);
151  }
152 
153  dynZones.append(zoneI);
154  dynSizes[zoneI]++;
155  }
156 
157  // skip empty groups
158  label nZone = 0;
159  forAll(dynSizes, zoneI)
160  {
161  if (dynSizes[zoneI])
162  {
163  if (nZone != zoneI)
164  {
165  dynSizes[nZone] = dynSizes[zoneI];
166  }
167  nZone++;
168  }
169  }
170  // truncate addressed size
171  dynSizes.setCapacity(nZone);
172 
173  // transfer to normal lists
174  points_.transfer(dynPoints);
175  zoneIds_.transfer(dynZones);
176  sizes_.transfer(dynSizes);
177 
178  return true;
179 }
180 
181 
182 // ************************************************************************* //
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
clear
UEqn clear()
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::IFstream
Input from file stream.
Definition: IFstream.H:81
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
IStringStream.H
Map.H
TRIsurfaceFormat.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:59
Foam::DynamicList::setCapacity
void setCapacity(const label)
Alter the size of the underlying storage.
Definition: DynamicListI.H:109
IOmanip.H
Istream and Ostream manipulators taking arguments.
IFstream.H
Foam::FatalError
error FatalError
Foam::fileFormats::TRIsurfaceFormatCore::read
bool read(const fileName &)
Definition: TRIsurfaceFormatCore.C:57
Foam::IStringStream
Input from memory buffer stream.
Definition: IStringStream.H:49
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::HashTable< label >
Foam::fileFormats::TRIsurfaceFormatCore::~TRIsurfaceFormatCore
~TRIsurfaceFormatCore()
Destructor.
Definition: TRIsurfaceFormatCore.C:50
Foam::fileFormats::TRIsurfaceFormatCore::TRIsurfaceFormatCore
TRIsurfaceFormatCore(const TRIsurfaceFormatCore &)
Disallow default bitwise copy construct.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::Vector< scalar >
Foam::line
A line primitive.
Definition: line.H:56
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::zero
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47