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-2022 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 void Foam::OBJstream::writeAndCheck(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 
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
65 (
66  const fileName& pathname,
67  IOstreamOption streamOpt
68 )
69 :
70  OFstream(pathname, streamOpt),
71  startOfLine_(true),
72  nVertices_(0)
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 {
80  writeAndCheck(c);
81  return *this;
82 }
83 
84 
85 Foam::Ostream& Foam::OBJstream::write(const char* str)
86 {
87  for (const char* iter = str; *iter; ++iter)
88  {
89  writeAndCheck(*iter);
90  }
91  return *this;
92 }
93 
94 
96 {
97  return writeQuoted(str, false);
98 }
99 
100 
102 {
103  return writeQuoted(str, true);
104 }
105 
106 
108 (
109  const std::string& str,
110  const bool quoted
111 )
112 {
113  if (!quoted)
114  {
115  // Output unquoted, only advance line number on newline
116  for (auto iter = str.cbegin(); iter != str.cend(); ++iter)
117  {
118  writeAndCheck(*iter);
119  }
120  return *this;
121  }
122 
123  // Output with surrounding quotes and backslash escaping
124  OFstream::write(static_cast<char>(token::DQUOTE));
125 
126  unsigned backslash = 0;
127  for (auto iter = str.cbegin(); iter != str.cend(); ++iter)
128  {
129  const char c = *iter;
130 
131  if (c == '\\')
132  {
133  ++backslash;
134  continue; // only output after escaped character is known
135  }
136  else if (c == token::NL)
137  {
138  ++lineNumber_;
139  ++backslash; // backslash escape for newline
140  }
141  else if (c == token::DQUOTE)
142  {
143  ++backslash; // backslash escape for quote
144  }
145 
146  // output all pending backslashes
147  while (backslash)
148  {
149  OFstream::write('\\');
150  --backslash;
151  }
152 
153  writeAndCheck(c);
154  }
155 
156  // silently drop any trailing backslashes
157  // they would otherwise appear like an escaped end-quote
158  OFstream::write(static_cast<char>(token::DQUOTE));
159 
160  return *this;
161 }
162 
163 
165 {
166  write('v') << ' ' << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
167  return *this;
168 }
169 
170 
172 {
173  write(p);
174  OFstream::write("vn ") << n.x() << ' ' << n.y() << ' ' << n.z() << nl;
175  return *this;
176 }
177 
178 
180 {
181  for (const point& p : points)
182  {
183  write('v') << ' ' << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
184  }
185  return *this;
186 }
187 
188 
189 Foam::Ostream& Foam::OBJstream::write(const edge& e, const UList<point>& points)
190 {
191  write(points[e.first()]);
192  write(points[e.second()]);
193  write('l') << ' ' << nVertices_-1 << ' ' << nVertices_ << nl;
194  return *this;
195 }
196 
197 
199 {
200  write(ln.first());
201  write(ln.second());
202  write('l') << ' ' << nVertices_-1 << ' ' << nVertices_ << nl;
203  return *this;
204 }
205 
206 
208 (
209  const linePointRef& ln,
210  const vector& n0,
211  const vector& n1
212 )
213 {
214  write(ln.first(), n0);
215  write(ln.second(), n1);
216  write('l') << ' ' << nVertices_-1 << ' ' << nVertices_ << nl;
217  return *this;
218 }
219 
220 
222 (
223  const point& p0,
224  const point& p1
225 )
226 {
227  write(p0);
228  write(p1);
229  write('l') << ' ' << nVertices_-1 << ' ' << nVertices_ << nl;
230  return *this;
231 }
232 
233 
235 (
236  const triPointRef& f,
237  const bool lines
238 )
239 {
240  const label start = nVertices_+1; // 1-offset for obj included here
241  write(f.a());
242  write(f.b());
243  write(f.c());
244  if (lines)
245  {
246  write('l');
247  for (int i = 0; i < 3; ++i)
248  {
249  write(' ') << i+start;
250  }
251  write(' ') << start << '\n';
252  }
253  else
254  {
255  write('f');
256  for (int i = 0; i < 3; ++i)
257  {
258  write(' ') << i+start;
259  }
260  write('\n');
261  }
262  return *this;
263 }
264 
265 
267 (
268  const UList<point>& points,
269  const bool lines
270 )
271 {
272  const label start = nVertices_+1; // 1-offset for obj included here
273 
274  write(points);
275 
276  if (lines)
277  {
278  write('l');
279  forAll(points, i)
280  {
281  write(' ') << i+start;
282  }
283  write(' ') << start << '\n';
284  }
285  else
286  {
287  write('f');
288  forAll(points, i)
289  {
290  write(' ') << i+start;
291  }
292  write('\n');
293  }
294  return *this;
295 }
296 
297 
299 (
300  const face& f,
301  const UList<point>& points,
302  const bool lines
303 )
304 {
305  const label start = nVertices_+1; // 1-offset for obj included here
306 
307  for (const label fp : f)
308  {
309  write(points[fp]);
310  }
311  if (lines)
312  {
313  write('l');
314  forAll(f, i)
315  {
316  write(' ') << i+start;
317  }
318  write(' ') << start << '\n';
319  }
320  else
321  {
322  write('f');
323  forAll(f, i)
324  {
325  write(' ') << i+start;
326  }
327  write('\n');
328  }
329  return *this;
330 }
331 
332 
334 (
335  const UList<face>& faces,
336  const pointField& points,
337  const bool lines
338 )
339 {
340  primitivePatch pp(SubList<face>(faces), points);
341 
342  const label start = nVertices_+1; // 1-offset for obj included here
343 
344  write(pp.localPoints());
345 
346  if (lines)
347  {
348  for (const edge& e : pp.edges())
349  {
350  write('l') << ' '
351  << e.first()+start << ' '
352  << e.second()+start << nl;
353  }
354  }
355  else
356  {
357  for (const face& f : pp.localFaces())
358  {
359  write('f');
360  for (const label fp : f)
361  {
362  write(' ') << fp+start;
363  }
364  write('\n');
365  }
366  }
367  return *this;
368 }
369 
370 
372 (
373  const UList<edge>& edges,
374  const UList<point>& points,
375  const bool compact
376 )
377 {
378  if (compact)
379  {
380  // Code similar to PrimitivePatch::calcMeshData()
381  // Unsorted version
382 
383  label objPointId = nVertices_+1; // 1-offset for obj included here
384 
385  Map<label> markedPoints(2*edges.size());
386 
387  for (const edge& e : edges)
388  {
389  if (markedPoints.insert(e.first(), objPointId))
390  {
391  write(points[e.first()]);
392  ++objPointId;
393  }
394  if (markedPoints.insert(e.second(), objPointId))
395  {
396  write(points[e.second()]);
397  ++objPointId;
398  }
399  }
400 
401  for (const edge& e : edges)
402  {
403  write('l') << ' '
404  << markedPoints[e.first()] << ' '
405  << markedPoints[e.second()] << nl;
406  }
407  }
408  else
409  {
410  const label start = nVertices_+1; // 1-offset for obj included here
411 
412  write(points);
413 
414  for (const edge& e : edges)
415  {
416  write('l') << ' '
417  << e.first()+start << ' '
418  << e.second()+start << nl;
419  }
420  }
421 
422  return *this;
423 }
424 
425 
427 (
428  const treeBoundBox& bb,
429  const bool lines
430 )
431 {
432  const label start = nVertices_+1; // 1-offset for obj included here
433 
434  write(bb.points());
435 
436  if (lines)
437  {
438  for (const edge& e : treeBoundBox::edges)
439  {
440  write('l') << ' '
441  << e.first()+start << ' '
442  << e.second()+start << nl;
443  }
444  }
445  else
446  {
447  for (const face& f : treeBoundBox::faces)
448  {
449  write('f');
450  for (const label fp : f)
451  {
452  write(' ') << fp+start;
453  }
454  write('\n');
455  }
456  }
457 
458  return *this;
459 }
460 
461 
462 // ************************************************************************* //
Ostream & writeFace(const UList< point > &points, const bool lines=true)
Write face loop points with lines/filled-polygon.
Definition: OBJstream.C:260
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:71
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
A class for handling file names.
Definition: fileName.H:71
const Field< point_type > & localPoints() const
Return pointField of points in patch.
OBJstream(const fileName &pathname, IOstreamOption streamOpt=IOstreamOption())
Construct from pathname (ASCII, uncompressed)
Definition: OBJstream.C:58
Double quote.
Definition: token.H:138
Output to file stream, using an OSstream.
Definition: OFstream.H:49
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
Newline [isspace].
Definition: token.H:127
A simple container for options an IOstream can normally have.
static const faceList faces
Face to point addressing, using octant corner points.
Definition: treeBoundBox.H:174
virtual Ostream & writeQuoted(const std::string &str, const bool quoted=true)
Write std::string surrounded by quotes.
Definition: OBJstream.C:101
os writeQuoted(("# "+outputName+"\), false)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
defineTypeName(manifoldCellsMeshObject)
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
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
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:99
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:55
Ostream & writeLine(const point &p0, const point &p1)
Write line joining two points.
Definition: OBJstream.C:215
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
virtual bool write(const token &tok)
Write token to stream or otherwise handle it.
Definition: OSstream.C:29
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.