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  Copyright (C) 2020-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "OBJstream.H"
30 #include "primitivePatch.H"
31 #include "treeBoundBox.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeName(OBJstream);
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 inline void Foam::OBJstream::vertex_state(const char c)
44 {
45  if (c == '\n')
46  {
47  startOfLine_ = true;
48  }
49  else if (startOfLine_)
50  {
51  startOfLine_ = false;
52  if (c == 'v')
53  {
54  ++nVertices_;
55  }
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const fileName& pathname,
65  IOstreamOption streamOpt
66 )
67 :
68  OFstream(pathname, streamOpt),
69  startOfLine_(true),
70  nVertices_(0)
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
79  vertex_state(c);
80  return *this;
81 }
82 
83 
85 (
86  const char* str,
87  std::streamsize len,
88  const bool quoted
89 )
90 {
91  OFstream::writeQuoted(str, len, quoted);
92 
93  // NOTE:
94  // Since vertex_state() handling only tracks newline followed by
95  // an initial 'v' (vertex) character, it probably should not be possible
96  // to triggered from within a quoted string.
97 
98  if (str && len >= 0)
99  {
100  if (quoted) vertex_state(0); // Begin quote: reset state
101 
102  const char* last = (str + len);
103 
104  // std::for_each
105  for (const char* iter = str; iter != last; ++iter)
106  {
107  vertex_state(*iter);
108  }
109 
110  if (quoted) vertex_state(0); // End quote
111  }
112 
113  return *this;
114 }
115 
116 
117 Foam::Ostream& Foam::OBJstream::write(const char* str)
118 {
119  OFstream::write(str);
120  for (const char* iter = str; *iter; ++iter)
121  {
122  vertex_state(*iter);
123  }
124  return *this;
125 }
126 
129 {
130  return writeQuoted(str.data(), str.size(), false);
131 }
132 
134 Foam::Ostream& Foam::OBJstream::write(const std::string& str)
135 {
136  return writeQuoted(str.data(), str.size(), true);
137 }
138 
139 
140 Foam::Ostream& Foam::OBJstream::writeComment(const std::string& str)
141 {
142  if (!startOfLine_)
143  {
144  OFstream::write('\n');
145  startOfLine_ = true;
146  }
147 
148  if (str.empty())
149  {
150  OFstream::write("#\n");
151  startOfLine_ = true;
152  return *this;
153  }
154 
155  const char* iter = str.data();
156  const char* last = (str.data() + str.size());
157 
158  // std::for_each
159  for (; iter != last; ++iter)
160  {
161  const char c = *iter;
162 
163  if (startOfLine_)
164  {
165  OFstream::write("# ");
166  startOfLine_ = false;
167  }
168 
169  startOfLine_ = (c == '\n');
171  }
172 
173  if (!startOfLine_)
174  {
175  OFstream::write('\n');
176  startOfLine_ = true;
177  }
178 
179  return *this;
180 }
181 
182 
184 {
185  write('v') << ' ' << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
186  return *this;
187 }
188 
189 
191 {
192  write(p);
193  OFstream::write("vn ") << n.x() << ' ' << n.y() << ' ' << n.z() << nl;
194  return *this;
195 }
196 
197 
199 {
200  for (const point& p : points)
201  {
202  write('v') << ' ' << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
203  }
204  return *this;
205 }
206 
207 
208 Foam::Ostream& Foam::OBJstream::write(const edge& e, const UList<point>& points)
209 {
210  write(points[e.first()]);
211  write(points[e.second()]);
212  write('l') << ' ' << nVertices_-1 << ' ' << nVertices_ << nl;
213  return *this;
214 }
215 
216 
218 {
219  write(ln.first());
220  write(ln.second());
221  write('l') << ' ' << nVertices_-1 << ' ' << nVertices_ << nl;
222  return *this;
223 }
224 
225 
227 (
228  const linePointRef& ln,
229  const vector& n0,
230  const vector& n1
231 )
232 {
233  write(ln.first(), n0);
234  write(ln.second(), n1);
235  write('l') << ' ' << nVertices_-1 << ' ' << nVertices_ << nl;
236  return *this;
237 }
238 
239 
241 (
242  const point& p0,
243  const point& p1
244 )
245 {
246  write(p0);
247  write(p1);
248  write('l') << ' ' << nVertices_-1 << ' ' << nVertices_ << nl;
249  return *this;
250 }
251 
252 
254 (
255  const triPointRef& f,
256  const bool lines
257 )
258 {
259  const label start = nVertices_+1; // 1-offset for obj included here
260  write(f.a());
261  write(f.b());
262  write(f.c());
263  if (lines)
264  {
265  write('l');
266  for (int i = 0; i < 3; ++i)
267  {
268  write(' ') << i+start;
269  }
270  write(' ') << start << '\n';
271  }
272  else
273  {
274  write('f');
275  for (int i = 0; i < 3; ++i)
276  {
277  write(' ') << i+start;
278  }
279  write('\n');
280  }
281  return *this;
282 }
283 
284 
286 (
287  const UList<point>& points,
288  const bool lines
289 )
290 {
291  const label start = nVertices_+1; // 1-offset for obj included here
292 
293  write(points);
294 
295  if (lines)
296  {
297  write('l');
298  forAll(points, i)
299  {
300  write(' ') << i+start;
301  }
302  write(' ') << start << '\n';
303  }
304  else
305  {
306  write('f');
307  forAll(points, i)
308  {
309  write(' ') << i+start;
310  }
311  write('\n');
312  }
313  return *this;
314 }
315 
316 
318 (
319  const face& f,
320  const UList<point>& points,
321  const bool lines
322 )
323 {
324  const label start = nVertices_+1; // 1-offset for obj included here
325 
326  for (const label fp : f)
327  {
328  write(points[fp]);
329  }
330  if (lines)
331  {
332  write('l');
333  forAll(f, i)
334  {
335  write(' ') << i+start;
336  }
337  write(' ') << start << '\n';
338  }
339  else
340  {
341  write('f');
342  forAll(f, i)
343  {
344  write(' ') << i+start;
345  }
346  write('\n');
347  }
348  return *this;
349 }
350 
351 
353 (
354  const UList<face>& faces,
355  const pointField& points,
356  const bool lines
357 )
358 {
359  primitivePatch pp(SubList<face>(faces), points);
360 
361  const label start = nVertices_+1; // 1-offset for obj included here
362 
363  write(pp.localPoints());
364 
365  if (lines)
366  {
367  for (const edge& e : pp.edges())
368  {
369  write('l') << ' '
370  << e.first()+start << ' '
371  << e.second()+start << nl;
372  }
373  }
374  else
375  {
376  for (const face& f : pp.localFaces())
377  {
378  write('f');
379  for (const label fp : f)
380  {
381  write(' ') << fp+start;
382  }
383  write('\n');
384  }
385  }
386  return *this;
387 }
388 
389 
391 (
392  const UList<edge>& edges,
393  const UList<point>& points,
394  const bool compact
395 )
396 {
397  if (compact)
398  {
399  // Code similar to PrimitivePatch::calcMeshData()
400  // Unsorted version
401 
402  label objPointId = nVertices_+1; // 1-offset for obj included here
403 
404  Map<label> markedPoints(2*edges.size());
405 
406  for (const edge& e : edges)
407  {
408  if (markedPoints.insert(e.first(), objPointId))
409  {
410  write(points[e.first()]);
411  ++objPointId;
412  }
413  if (markedPoints.insert(e.second(), objPointId))
414  {
415  write(points[e.second()]);
416  ++objPointId;
417  }
418  }
419 
420  for (const edge& e : edges)
421  {
422  write('l') << ' '
423  << markedPoints[e.first()] << ' '
424  << markedPoints[e.second()] << nl;
425  }
426  }
427  else
428  {
429  const label start = nVertices_+1; // 1-offset for obj included here
430 
431  write(points);
432 
433  for (const edge& e : edges)
434  {
435  write('l') << ' '
436  << e.first()+start << ' '
437  << e.second()+start << nl;
438  }
439  }
440 
441  return *this;
442 }
443 
444 
446 (
447  const treeBoundBox& bb,
448  const bool lines
449 )
450 {
451  const label start = nVertices_+1; // 1-offset for obj included here
452 
453  write(bb.points());
454 
455  if (lines)
456  {
457  for (const edge& e : treeBoundBox::edges)
458  {
459  write('l') << ' '
460  << e.first()+start << ' '
461  << e.second()+start << nl;
462  }
463  }
464  else
465  {
466  for (const face& f : treeBoundBox::faces)
467  {
468  write('f');
469  for (const label fp : f)
470  {
471  write(' ') << fp+start;
472  }
473  write('\n');
474  }
475  }
476 
477  return *this;
478 }
479 
480 
481 // ************************************************************************* //
Ostream & writeFace(const UList< point > &points, const bool lines=true)
Write face loop points with lines/filled-polygon.
Definition: OBJstream.C:279
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:217
A line primitive.
Definition: line.H:52
const Field< point_type > & localPoints() const
Return pointField of points in patch.
virtual Ostream & write(const char c) override
Write character.
Definition: OBJstream.C:69
virtual Ostream & writeQuoted(const char *str, std::streamsize len, const bool quoted=true) override
Write character/string content, with/without surrounding quotes.
Definition: OSstream.C:101
OBJstream(const fileName &pathname, IOstreamOption streamOpt=IOstreamOption())
Construct from pathname (ASCII, uncompressed)
Definition: OBJstream.C:56
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static const faceList faces
Face to point addressing, using octant corner points.
Definition: treeBoundBox.H:174
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
defineTypeName(manifoldCellsMeshObject)
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Ostream & writeComment(const std::string &str)
Write comment (with &#39;# &#39; prefix)
Definition: OBJstream.C:133
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const pointField & points
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual Ostream & writeQuoted(const char *str, std::streamsize len, const bool quoted=true) override
Write character/string content, with/without surrounding quotes.
Definition: OBJstream.C:78
static const edgeList edges
Edge to point addressing, using octant corner points.
Definition: treeBoundBox.H:179
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
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:105
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:1237
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
virtual bool write(const token &tok) override
Write token to stream or otherwise handle it.
Definition: OSstream.C:29
Ostream & writeLine(const point &p0, const point &p1)
Write line joining two points.
Definition: OBJstream.C:234
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
labelList f(nPoints)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
const dimensionedScalar c
Speed of light in a vacuum.
label n
volScalarField & p
const volScalarField & p0
Definition: EEqn.H:36
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.