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-------------------------------------------------------------------------------
11License
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
35namespace Foam
36{
38}
39
40
41// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42
43void 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
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 ") << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
167 return *this;
168}
169
170
172{
173 write(pt);
174 OFstream::write("vn ") << n.x() << ' ' << n.y() << ' ' << n.z() << nl;
175 return *this;
176}
177
178
180{
181 write(points[e[0]]);
182 write(points[e[1]]);
183 write("l ") << nVertices_-1 << ' ' << nVertices_ << nl;
184 return *this;
185}
186
187
189{
190 write(ln.start());
191 write(ln.end());
192 write("l ") << nVertices_-1 << ' ' << nVertices_ << nl;
193 return *this;
194}
195
196
198(
199 const linePointRef& ln,
200 const vector& n0,
201 const vector& n1
202)
203{
204 write(ln.start(), n0);
205 write(ln.end(), n1);
206 write("l ") << nVertices_-1 << ' ' << nVertices_ << nl;
207 return *this;
208}
209
210
212(
213 const triPointRef& f,
214 const bool lines
215)
216{
217 const label start = nVertices_+1; // 1-offset for obj included here
218 write(f.a());
219 write(f.b());
220 write(f.c());
221 if (lines)
222 {
223 write('l');
224 for (int i = 0; i < 3; i++)
225 {
226 write(' ') << i+start;
227 }
228 write(' ') << start << '\n';
229 }
230 else
231 {
232 write('f');
233 for (int i = 0; i < 3; i++)
234 {
235 write(' ') << i+start;
236 }
237 write('\n');
238 }
239 return *this;
240}
241
242
244(
245 const face& f,
246 const UList<point>& points,
247 const bool lines
248)
249{
250 const label start = nVertices_+1; // 1-offset for obj included here
251 forAll(f, i)
252 {
253 write(points[f[i]]);
254 }
255 if (lines)
256 {
257 write('l');
258 forAll(f, i)
259 {
260 write(' ') << i+start;
261 }
262 write(' ') << start << '\n';
263 }
264 else
265 {
266 write('f');
267 forAll(f, i)
268 {
269 write(' ') << i+start;
270 }
271 write('\n');
272 }
273 return *this;
274}
275
276
278(
279 const UList<face>& faces,
280 const pointField& points,
281 const bool lines
282)
283{
285
286 const pointField& localPoints = pp.localPoints();
287 const faceList& localFaces = pp.localFaces();
288
289 const label start = nVertices_+1; // 1-offset for obj included here
290
291 forAll(localPoints, i)
292 {
293 write(localPoints[i]);
294 }
295
296 if (lines)
297 {
298 const edgeList& edges = pp.edges();
299 forAll(edges, edgeI)
300 {
301 const edge& e = edges[edgeI];
302
303 write("l ") << e[0]+start << ' ' << e[1]+start << nl;
304 }
305 }
306 else
307 {
308 forAll(localFaces, facei)
309 {
310 const face& f = localFaces[facei];
311 write('f');
312 forAll(f, i)
313 {
314 write(' ') << f[i]+start;
315 }
316 write('\n');
317 }
318 }
319 return *this;
320}
321
322
324(
325 const UList<edge>& edges,
326 const UList<point>& points,
327 const bool compact
328)
329{
330 if (compact)
331 {
332 // Code similar to PrimitivePatch::calcMeshData()
333 // Unsorted version
334
335 label objPointId = nVertices_+1; // 1-offset for obj included here
336
337 Map<label> markedPoints(2*edges.size());
338 forAll(edges, edgei)
339 {
340 const edge& e = edges[edgei];
341
342 if (markedPoints.insert(e[0], objPointId))
343 {
344 write(points[e[0]]);
345 ++objPointId;
346 }
347 if (markedPoints.insert(e[1], objPointId))
348 {
349 write(points[e[1]]);
350 ++objPointId;
351 }
352 }
353
354 forAll(edges, edgei)
355 {
356 const edge& e = edges[edgei];
357
358 write("l ")
359 << markedPoints[e[0]] << ' '
360 << markedPoints[e[1]] << nl;
361 }
362 }
363 else
364 {
365 const label start = nVertices_+1; // 1-offset for obj included here
366
367 forAll(points, i)
368 {
369 write(points[i]);
370 }
371
372 forAll(edges, edgei)
373 {
374 const edge& e = edges[edgei];
375
376 write("l ")
377 << e[0]+start << ' ' << e[1]+start << nl;
378 }
379 }
380
381 return *this;
382}
383
384
386(
387 const treeBoundBox& bb,
388 const bool lines
389)
390{
391 const label start = nVertices_+1; // 1-offset for obj included here
392
394 forAll(points, i)
395 {
396 write(points[i]);
397 }
398
399 if (lines)
400 {
402 {
403 const edge& e = treeBoundBox::edges[edgei];
404
405 write("l ") << e[0]+start << ' ' << e[1]+start << nl;
406 }
407 }
408 else
409 {
411 {
412 const face& f = treeBoundBox::faces[facei];
413
414 write('f');
415 forAll(f, i)
416 {
417 write(' ') << f[i]+start;
418 }
419 write('\n');
420 }
421 }
422
423 return *this;
424}
425
426
427// ************************************************************************* //
label n
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
The IOstreamOption is a simple container for options an IOstream can normally have.
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
OFstream that keeps track of vertices.
Definition: OBJstream.H:61
virtual Ostream & writeQuoted(const std::string &str, const bool quoted=true)
Write std::string surrounded by quotes.
Definition: OBJstream.C:108
Output to file stream, using an OSstream.
Definition: OFstream.H:57
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of faces which address into the list of points.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
A List obtained as a section of another List.
Definition: SubList.H:70
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A class for handling file names.
Definition: fileName.H:76
virtual bool write()
Write the output fields.
A line primitive.
Definition: line.H:68
@ DQUOTE
Double quote.
Definition: token.H:135
@ NL
Newline [isspace].
Definition: token.H:124
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:154
static const faceList faces
Face to point addressing.
Definition: treeBoundBox.H:151
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:92
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:80
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
os writeQuoted(("# "+outputName+"\n"), false)
const pointField & points
Namespace for OpenFOAM.
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: MSwindows.C:933
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
runTime write()
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333