OBJstream.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) 2012-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "OBJstream.H"
29 #include "primitivePatch.H"
30 #include "treeBoundBox.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 defineTypeNameAndDebug(OBJstream, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 void Foam::OBJstream::writeAndCheck(const char c)
43 {
44  if (c == '\n')
45  {
46  startOfLine_ = true;
47  }
48  else if (startOfLine_)
49  {
50  startOfLine_ = false;
51  if (c == 'v')
52  {
53  ++nVertices_;
54  }
55  }
56 
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
64 (
65  const fileName& pathname,
68  compressionType compression
69 )
70 :
71  OFstream(pathname, format, version, compression),
72  startOfLine_(true),
73  nVertices_(0)
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  writeAndCheck(c);
82  return *this;
83 }
84 
85 
87 {
88  for (const char* iter = str; *iter; ++iter)
89  {
90  writeAndCheck(*iter);
91  }
92  return *this;
93 }
94 
95 
97 {
98  return writeQuoted(str, false);
99 }
100 
101 
103 {
104  return writeQuoted(str, true);
105 }
106 
107 
109 (
110  const std::string& str,
111  const bool quoted
112 )
113 {
114  if (!quoted)
115  {
116  // Output unquoted, only advance line number on newline
117  for (auto iter = str.cbegin(); iter != str.cend(); ++iter)
118  {
119  writeAndCheck(*iter);
120  }
121  return *this;
122  }
123 
124 
125  OFstream::write(static_cast<char>(token::BEGIN_STRING));
126 
127  unsigned backslash = 0;
128  for (auto iter = str.cbegin(); iter != str.cend(); ++iter)
129  {
130  const char c = *iter;
131 
132  if (c == '\\')
133  {
134  ++backslash;
135  continue; // only output after escaped character is known
136  }
137  else if (c == token::NL)
138  {
139  ++lineNumber_;
140  ++backslash; // backslash escape for newline
141  }
142  else if (c == token::END_STRING)
143  {
144  ++backslash; // backslash escape for quote
145  }
146 
147  // output all pending backslashes
148  while (backslash)
149  {
150  OFstream::write('\\');
151  --backslash;
152  }
153 
154  writeAndCheck(c);
155  }
156 
157  // silently drop any trailing backslashes
158  // they would otherwise appear like an escaped end-quote
159  OFstream::write(static_cast<char>(token::END_STRING));
160 
161  return *this;
162 }
163 
164 
166 {
167  write("v ") << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
168  return *this;
169 }
170 
171 
173 {
174  write(pt);
175  OFstream::write("vn ") << n.x() << ' ' << n.y() << ' ' << n.z() << nl;
176  return *this;
177 }
178 
179 
181 {
182  write(points[e[0]]);
183  write(points[e[1]]);
184  write("l ") << nVertices_-1 << ' ' << nVertices_ << nl;
185  return *this;
186 }
187 
188 
190 {
191  write(ln.start());
192  write(ln.end());
193  write("l ") << nVertices_-1 << ' ' << nVertices_ << nl;
194  return *this;
195 }
196 
197 
199 (
200  const linePointRef& ln,
201  const vector& n0,
202  const vector& n1
203 )
204 {
205  write(ln.start(), n0);
206  write(ln.end(), n1);
207  write("l ") << nVertices_-1 << ' ' << nVertices_ << nl;
208  return *this;
209 }
210 
211 
213 (
214  const triPointRef& f,
215  const bool lines
216 )
217 {
218  const label start = nVertices_+1; // 1-offset for obj included here
219  write(f.a());
220  write(f.b());
221  write(f.c());
222  if (lines)
223  {
224  write('l');
225  for (int i = 0; i < 3; i++)
226  {
227  write(' ') << i+start;
228  }
229  write(' ') << start << '\n';
230  }
231  else
232  {
233  write('f');
234  for (int i = 0; i < 3; i++)
235  {
236  write(' ') << i+start;
237  }
238  write('\n');
239  }
240  return *this;
241 }
242 
243 
245 (
246  const face& f,
247  const UList<point>& points,
248  const bool lines
249 )
250 {
251  const label start = nVertices_+1; // 1-offset for obj included here
252  forAll(f, i)
253  {
254  write(points[f[i]]);
255  }
256  if (lines)
257  {
258  write('l');
259  forAll(f, i)
260  {
261  write(' ') << i+start;
262  }
263  write(' ') << start << '\n';
264  }
265  else
266  {
267  write('f');
268  forAll(f, i)
269  {
270  write(' ') << i+start;
271  }
272  write('\n');
273  }
274  return *this;
275 }
276 
277 
279 (
280  const UList<face>& faces,
281  const pointField& points,
282  const bool lines
283 )
284 {
285  SubList<face> allFcs(faces, faces.size());
286 
287  primitivePatch pp(allFcs, points);
288 
289  const pointField& localPoints = pp.localPoints();
290  const faceList& localFaces = pp.localFaces();
291 
292  const label start = nVertices_+1; // 1-offset for obj included here
293 
294  forAll(localPoints, i)
295  {
296  write(localPoints[i]);
297  }
298 
299  if (lines)
300  {
301  const edgeList& edges = pp.edges();
302  forAll(edges, edgeI)
303  {
304  const edge& e = edges[edgeI];
305 
306  write("l ") << e[0]+start << ' ' << e[1]+start << nl;
307  }
308  }
309  else
310  {
311  forAll(localFaces, facei)
312  {
313  const face& f = localFaces[facei];
314  write('f');
315  forAll(f, i)
316  {
317  write(' ') << f[i]+start;
318  }
319  write('\n');
320  }
321  }
322  return *this;
323 }
324 
325 
327 (
328  const UList<edge>& edges,
329  const UList<point>& points,
330  const bool compact
331 )
332 {
333  if (compact)
334  {
335  // Code similar to PrimitivePatch::calcMeshData()
336  // Unsorted version
337 
338  label objPointId = nVertices_+1; // 1-offset for obj included here
339 
340  Map<label> markedPoints(2*edges.size());
341  forAll(edges, edgei)
342  {
343  const edge& e = edges[edgei];
344 
345  if (markedPoints.insert(e[0], objPointId))
346  {
347  write(points[e[0]]);
348  ++objPointId;
349  }
350  if (markedPoints.insert(e[1], objPointId))
351  {
352  write(points[e[1]]);
353  ++objPointId;
354  }
355  }
356 
357  forAll(edges, edgei)
358  {
359  const edge& e = edges[edgei];
360 
361  write("l ")
362  << markedPoints[e[0]] << ' '
363  << markedPoints[e[1]] << nl;
364  }
365  }
366  else
367  {
368  const label start = nVertices_+1; // 1-offset for obj included here
369 
370  forAll(points, i)
371  {
372  write(points[i]);
373  }
374 
375  forAll(edges, edgei)
376  {
377  const edge& e = edges[edgei];
378 
379  write("l ")
380  << e[0]+start << ' ' << e[1]+start << nl;
381  }
382  }
383 
384  return *this;
385 }
386 
387 
389 (
390  const treeBoundBox& bb,
391  const bool lines
392 )
393 {
394  const label start = nVertices_+1; // 1-offset for obj included here
395 
396  pointField points(bb.points());
397  forAll(points, i)
398  {
399  write(points[i]);
400  }
401 
402  if (lines)
403  {
405  {
406  const edge& e = treeBoundBox::edges[edgei];
407 
408  write("l ") << e[0]+start << ' ' << e[1]+start << nl;
409  }
410  }
411  else
412  {
414  {
415  const face& f = treeBoundBox::faces[facei];
416 
417  write('f');
418  forAll(f, i)
419  {
420  write(' ') << f[i]+start;
421  }
422  write('\n');
423  }
424  }
425 
426  return *this;
427 }
428 
429 
430 // ************************************************************************* //
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:78
Foam::OSstream::write
virtual bool write(const token &tok)
Write token to stream or otherwise handle it.
Definition: OSstream.C:36
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::treeBoundBox::points
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:92
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::OBJstream::writeQuoted
virtual Ostream & writeQuoted(const std::string &str, const bool quoted=true)
Write std::string surrounded by quotes.
Definition: OBJstream.C:109
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:87
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::token::END_STRING
End string with double quote.
Definition: token.H:138
Foam::Map< label >
Foam::treeBoundBox::edges
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:153
Foam::OBJstream::write
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:79
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:90
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:62
n
label n
Definition: TABSMDCalcMethod2.H:31
format
word format(conversionProperties.get< word >("format"))
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:62
Foam::IOstreamOption::versionNumber
Representation of a major/minor version number.
Definition: IOstreamOption.H:79
Foam::Field< vector >
treeBoundBox.H
Foam::treeBoundBox::faces
static const faceList faces
Face to point addressing.
Definition: treeBoundBox.H:150
Foam::OBJstream::OBJstream
OBJstream(const fileName &pathname, streamFormat format=ASCII, versionNumber version=currentVersion, compressionType compression=UNCOMPRESSED)
Construct from pathname.
Definition: OBJstream.C:64
Foam::IOstreamOption::streamFormat
streamFormat
Data format (ascii | binary)
Definition: IOstreamOption.H:64
Foam::token::NL
Newline [isspace].
Definition: token.H:114
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::foamVersion::version
const std::string version
OpenFOAM version (name or stringified number) as a std::string.
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:99
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:84
Foam::nl
constexpr char nl
Definition: Ostream.H:372
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< face >
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::line
A line primitive.
Definition: line.H:59
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
Foam::IOstreamOption::compressionType
compressionType
Compression treatment (UNCOMPRESSED | COMPRESSED)
Definition: IOstreamOption.H:71
Foam::ln
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: MSwindows.C:915
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
primitivePatch.H
Foam::token::BEGIN_STRING
Begin string with double quote.
Definition: token.H:137
OBJstream.H
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:90