OBJedgeFormat.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "OBJedgeFormat.H"
27 #include "clock.H"
28 #include "IFstream.H"
29 #include "IStringStream.H"
30 #include "Ostream.H"
31 #include "OFstream.H"
32 #include "ListOps.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 Foam::fileFormats::OBJedgeFormat::OBJedgeFormat
37 (
38  const fileName& filename
39 )
40 {
41  read(filename);
42 }
43 
44 
45 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
46 
47 void Foam::fileFormats::OBJedgeFormat::readVertices
48 (
49  const string& line,
50  string::size_type& endNum,
51  DynamicList<label>& dynVertices
52 )
53 {
54  dynVertices.clear();
55  while (true)
56  {
57  string::size_type startNum =
58  line.find_first_not_of(' ', endNum);
59 
60  if (startNum == string::npos)
61  {
62  break;
63  }
64 
65  endNum = line.find(' ', startNum);
66 
67  string vertexSpec;
68  if (endNum != string::npos)
69  {
70  vertexSpec = line.substr(startNum, endNum-startNum);
71  }
72  else
73  {
74  vertexSpec = line.substr(startNum, line.size() - startNum);
75  }
76 
77  string::size_type slashPos = vertexSpec.find('/');
78 
79  label vertI = 0;
80  if (slashPos != string::npos)
81  {
82  IStringStream intStream(vertexSpec.substr(0, slashPos));
83 
84  intStream >> vertI;
85  }
86  else
87  {
88  IStringStream intStream(vertexSpec);
89 
90  intStream >> vertI;
91  }
92  dynVertices.append(vertI - 1);
93  }
94 }
95 
96 
98 {
99  clear();
100 
101  IFstream is(filename);
102  if (!is.good())
103  {
105  << "Cannot read file " << filename
106  << exit(FatalError);
107  }
108 
109  DynamicList<point> dynPoints;
110  DynamicList<edge> dynEdges;
111  DynamicList<label> dynUsedPoints;
112 
113  DynamicList<label> dynVertices;
114 
115  while (is.good())
116  {
117  string line = this->getLineNoComment(is);
118 
119  // handle continuations
120  if (line[line.size()-1] == '\\')
121  {
122  line.substr(0, line.size()-1);
123  line += this->getLineNoComment(is);
124  }
125 
126  // Read first word
127  IStringStream lineStream(line);
128  word cmd;
129  lineStream >> cmd;
130 
131  if (cmd == "v")
132  {
133  scalar x, y, z;
134  lineStream >> x >> y >> z;
135 
136  dynPoints.append(point(x, y, z));
137  dynUsedPoints.append(-1);
138  }
139  else if (cmd == "l")
140  {
141  // Assume 'l' is followed by space.
142  string::size_type endNum = 1;
143 
144  readVertices
145  (
146  line,
147  endNum,
148  dynVertices
149  );
150 
151 
152  for (label i = 1; i < dynVertices.size(); i++)
153  {
154  edge edgeRead(dynVertices[i-1], dynVertices[i]);
155 
156  dynUsedPoints[edgeRead[0]] = edgeRead[0];
157  dynUsedPoints[edgeRead[1]] = edgeRead[1];
158 
159  dynEdges.append(edgeRead);
160  }
161  }
162  else if (cmd == "f")
163  {
164  // Support for faces with 2 vertices
165 
166  // Assume 'f' is followed by space.
167  string::size_type endNum = 1;
168 
169  readVertices
170  (
171  line,
172  endNum,
173  dynVertices
174  );
175 
176  if (dynVertices.size() == 2)
177  {
178  for (label i = 1; i < dynVertices.size(); i++)
179  {
180  edge edgeRead(dynVertices[i-1], dynVertices[i]);
181 
182  dynUsedPoints[edgeRead[0]] = edgeRead[0];
183  dynUsedPoints[edgeRead[1]] = edgeRead[1];
184 
185  dynEdges.append(edgeRead);
186  }
187  }
188  }
189  }
190 
191  // cull unused points
192  label nUsed = 0;
193 
194  forAll(dynPoints, pointi)
195  {
196  if (dynUsedPoints[pointi] >= 0)
197  {
198  if (nUsed != pointi)
199  {
200  dynPoints[nUsed] = dynPoints[pointi];
201  dynUsedPoints[pointi] = nUsed; // new position
202  }
203  ++nUsed;
204  }
205  }
206 
207  dynPoints.setSize(nUsed);
208 
209  // transfer to normal lists
210  storedPoints().transfer(dynPoints);
211 
212  // renumber edge vertices
213  if (nUsed != dynUsedPoints.size())
214  {
215  forAll(dynEdges, edgeI)
216  {
217  edge& e = dynEdges[edgeI];
218 
219  e[0] = dynUsedPoints[e[0]];
220  e[1] = dynUsedPoints[e[1]];
221  }
222  }
223  storedEdges().transfer(dynEdges);
224 
225  return true;
226 }
227 
228 
230 (
231  const fileName& filename,
232  const edgeMesh& mesh
233 )
234 {
235  const pointField& pointLst = mesh.points();
236  const edgeList& edgeLst = mesh.edges();
237 
238  OFstream os(filename);
239  if (!os.good())
240  {
242  << "Cannot open file for writing " << filename
243  << exit(FatalError);
244  }
245 
246 
247  os << "# Wavefront OBJ file written " << clock::dateTime().c_str() << nl
248  << "o " << os.name().lessExt().name() << nl
249  << nl
250  << "# points : " << pointLst.size() << nl
251  << "# lines : " << edgeLst.size() << nl;
252 
253  os << nl
254  << "# <points count=\"" << pointLst.size() << "\">" << nl;
255 
256  // Write vertex coords
257  forAll(pointLst, ptI)
258  {
259  const point& p = pointLst[ptI];
260 
261  os << "v " << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
262  }
263 
264  os << "# </points>" << nl
265  << nl
266  << "# <edges count=\"" << edgeLst.size() << "\">" << endl;
267 
268  // Write line connectivity
269  forAll(edgeLst, edgeI)
270  {
271  const edge& e = edgeLst[edgeI];
272 
273  os << "l " << (e[0] + 1) << " " << (e[1] + 1) << nl;
274  }
275  os << "# </edges>" << endl;
276 }
277 
278 
279 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
edgeList & storedEdges()
Non-const access to the edges.
Definition: edgeMeshI.H:34
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Output to file stream.
Definition: OFstream.H:82
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
virtual bool read(const fileName &)
Read from file.
Definition: OBJedgeFormat.C:97
const Cmpt & z() const
Definition: VectorI.H:87
static void write(const fileName &, const edgeMesh &)
Write surface mesh components by proxy.
static string getLineNoComment(IFstream &)
Read non-comment line.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
pointField & storedPoints()
Non-const access to global points.
Definition: edgeMeshI.H:28
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
const Cmpt & y() const
Definition: VectorI.H:81
Various functions to operate on Lists.
scalar y
void setSize(const label)
Alter the addressed list size.
Definition: DynamicListI.H:163
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
A class for handling words, derived from string.
Definition: word.H:59
static string dateTime()
Return the current wall-clock date/time as a string.
Definition: clock.C:57
const edgeList & edges() const
Return edges.
Definition: edgeMeshI.H:59
word name() const
Return file name (part beyond last /)
Definition: fileName.C:179
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
const Cmpt & x() const
Definition: VectorI.H:75
static const char nl
Definition: Ostream.H:265
const pointField & points() const
Return points.
Definition: edgeMeshI.H:53
Points connected by edges.
Definition: edgeMesh.H:69
Input from file stream.
Definition: IFstream.H:81
virtual void clear()
Clear all storage.
Definition: edgeMesh.C:169
vector point
Point is a vector.
Definition: point.H:41
Input from memory buffer stream.
Definition: IStringStream.H:49
fileName lessExt() const
Return file name without extension (part before last .)
Definition: fileName.C:268
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
volScalarField & p
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:224
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342