readTRI.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | foam-extend: Open Source CFD
4  \\ / O peration | Version: 3.2
5  \\ / A nd | Web: http://www.foam-extend.org
6  \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9  This file is part of foam-extend.
10 
11  foam-extend is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation, either version 3 of the License, or (at your
14  option) any later version.
15 
16  foam-extend is distributed in the hope that it will be useful, but
17  WITHOUT ANY WARRANTY; without even the implied warranty of
18  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19  General Public License for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  TRI (triangle) file reader. Comes out of e.g. AC3D.
26  lines are 9 floats (3 points, each 3 floats) followed by hex colour.
27  Is converted into regions: regions numbered from 0 up, each colour is
28  region.
29  Most of reading/stitching taken from STL reader.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "triSurface.H"
34 #include "STLpoint.H"
35 #include "SLList.H"
36 #include "IFstream.H"
37 #include "readHexLabel.H"
38 #include "stringList.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
46 
47 bool triSurface::readTRI(const fileName& TRIfileName)
48 {
49  IFstream TRIfile(TRIfileName);
50 
51  if (!TRIfile.good())
52  {
53  FatalErrorIn("triSurface::readTRI(const fileName&)")
54  << "Cannot read file " << TRIfileName
55  << exit(FatalError);
56  }
57 
58  SLList<STLpoint> STLpoints;
59  SLList<label> STLlabels;
60  HashTable<label, string> STLsolidNames;
61 
62  // Max region number so far
63  label maxRegion = 0;
64 
65  while(TRIfile)
66  {
67  string line = getLineNoComment(TRIfile);
68 
69  if (line.empty())
70  {
71  break;
72  }
73 
74  IStringStream lineStream(line);
75 
76  STLpoint p
77  (
78  readScalar(lineStream),
79  readScalar(lineStream),
80  readScalar(lineStream)
81  );
82 
83  if (!lineStream) break;
84 
85  STLpoints.append(p);
86  STLpoints.append
87  (
88  STLpoint
89  (
90  readScalar(lineStream),
91  readScalar(lineStream),
92  readScalar(lineStream)
93  )
94  );
95  STLpoints.append
96  (
97  STLpoint
98  (
99  readScalar(lineStream),
100  readScalar(lineStream),
101  readScalar(lineStream)
102  )
103  );
104 
105  // Region/colour in .tri file starts with 0x. Skip.
106 
107  char zero;
108  lineStream >> zero;
109 
110  word rawSolidName(lineStream);
111 
112  word solidName("patch" + rawSolidName(1, rawSolidName.size()-1));
113 
114  label region = -1;
115 
117  STLsolidNames.find(solidName);
118 
119  if (findName != STLsolidNames.end())
120  {
121  region = findName();
122  }
123  else
124  {
125  Pout<< "Mapping triangle colour 0" << rawSolidName
126  << " to region " << maxRegion << " name " << solidName
127  << endl;
128 
129  region = maxRegion++;
130  STLsolidNames.insert(solidName, region);
131  }
132  STLlabels.append(region);
133  }
134 
135 
136  pointField rawPoints(STLpoints.size());
137 
138  label i = 0;
139  for
140  (
141  SLList<STLpoint>::iterator iter = STLpoints.begin();
142  iter != STLpoints.end();
143  ++iter
144  )
145  {
146  rawPoints[i++] = *iter;
147  }
148 
149  setSize(STLlabels.size());
150 
151  label pointI = 0;
152  SLList<label>::iterator iter = STLlabels.begin();
153  forAll (*this, i)
154  {
155  operator[](i)[0] = pointI++;
156  operator[](i)[1] = pointI++;
157  operator[](i)[2] = pointI++;
158  operator[](i).region() = *iter;
159  ++iter;
160  }
161 
162  stitchTriangles(rawPoints);
163 
164  // Convert solidNames into regionNames
165  stringList names(STLsolidNames.toc());
166 
167  patches_.setSize(names.size());
168 
169  forAll(names, nameI)
170  {
171  patches_[nameI].name() = names[nameI];
172  patches_[nameI].geometricType() = "empty";
173  }
174 
175  return true;
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace Foam
182 
183 // ************************************************************************* //
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::LList< SLListBase, T >::append
void append(const T &a)
Add at tail of list.
Definition: LList.H:166
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::HashTable::toc
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:201
Foam::SLList
Non-intrusive singly-linked list.
Definition: SLList.H:47
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::HashTable::const_iterator
An STL-conforming const_iterator.
Definition: HashTable.H:470
Foam::HashTable::insert
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
readHexLabel.H
Read a hex label from an input stream.
Foam::triSurface::readTRI
bool readTRI(const fileName &)
Definition: readTRI.C:47
Foam::LList< SLListBase, T >::end
const iterator & end()
Definition: LList.H:273
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::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::LList< SLListBase, T >::begin
iterator begin()
Definition: LList.H:268
IFstream.H
Foam::FatalError
error FatalError
Foam::IStringStream
Input from memory buffer stream.
Definition: IStringStream.H:49
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::HashTable::find
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
Foam::List< labelledTri >::setSize
void setSize(const label)
Reset size of List.
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::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::triSurface::patches_
geometricSurfacePatchList patches_
Patch information (face ordering nFaces/startFace only used.
Definition: triSurface.H:83
SLList.H
Foam::STLpoint
A vertex point representation for STL files.
Definition: STLpoint.H:48
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::triSurface::stitchTriangles
bool stitchTriangles(const pointField &rawPoints, const scalar tol=SMALL, const bool verbose=false)
Function to stitch the triangles by removing duplicate points.
Definition: stitchTriangles.C:38
Foam::line
A line primitive.
Definition: line.H:56
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
stringList.H
Foam::triSurface::getLineNoComment
static string getLineNoComment(IFstream &)
Read non-comment line.
Definition: triSurface.C:172
Foam::zero
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47