OBJedgeFormat.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2018 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 "OBJedgeFormat.H"
30#include "clock.H"
31#include "Fstream.H"
32#include "Ostream.H"
33#include "stringOps.H"
34
35// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39
40// Token list with one of the following:
41// f v1 v2 v3 ...
42// f v1/vt1 v2/vt2 v3/vt3 ...
43// l v1 v2 v3 ...
44// l v1/vt1 v2/vt2 v3/vt3 ...
45static label readObjVertices
46(
47 const SubStrings<string>& tokens,
49)
50{
51 verts.clear();
52
53 bool first = true;
54 for (const auto& tok : tokens)
55 {
56 if (first)
57 {
58 // skip initial "f" or "l"
59 first = false;
60 continue;
61 }
62
63 const string vrtSpec(tok);
64 const auto slash = vrtSpec.find('/');
65
66 const label vertId =
67 (
68 slash != string::npos
69 ? readLabel(vrtSpec.substr(0, slash))
70 : readLabel(vrtSpec)
71 );
72
73 verts.append(vertId - 1);
74 }
75
76 return verts.size();
77}
78
79} // End namespace Foam
80
81
82
83// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84
86(
87 const fileName& filename
88)
89{
90 read(filename);
91}
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
97{
98 clear();
99
100 IFstream is(filename);
101 if (!is.good())
102 {
104 << "Cannot read file " << filename
105 << exit(FatalError);
106 }
107
108 DynamicList<point> dynPoints;
109 DynamicList<label> dynVerts;
110 DynamicList<edge> dynEdges;
111 DynamicList<label> dynUsedPoints;
112
113 while (is.good())
114 {
115 string line = this->getLineNoComment(is);
116
117 // Line continuations
118 while (line.removeEnd('\\'))
119 {
120 line += this->getLineNoComment(is);
121 }
122
124
125 // Require command and some arguments
126 if (tokens.size() < 2)
127 {
128 continue;
129 }
130
131 const word cmd = word::validate(tokens[0]);
132
133 if (cmd == "v")
134 {
135 // Vertex
136 // v x y z
137
138 dynPoints.append
139 (
140 point
141 (
142 readScalar(tokens[1]),
143 readScalar(tokens[2]),
144 readScalar(tokens[3])
145 )
146 );
147
148 dynUsedPoints.append(-1);
149 }
150 else if (cmd == "l")
151 {
152 // Line
153 // l v1 v2 v3 ...
154 // OR
155 // l v1/vt1 v2/vt2 v3/vt3 ...
156
157 readObjVertices(tokens, dynVerts);
158
159 for (label i = 1; i < dynVerts.size(); i++)
160 {
161 const edge e(dynVerts[i-1], dynVerts[i]);
162 dynEdges.append(e);
163
164 dynUsedPoints[e[0]] = e[0];
165 dynUsedPoints[e[1]] = e[1];
166 }
167 }
168 else if (cmd == "f")
169 {
170 // Face - support for faces with 2 vertices
171
172 // f v1 v2 v3 ...
173 // OR
174 // f v1/vt1 v2/vt2 v3/vt3 ...
175
176 if (readObjVertices(tokens, dynVerts) == 2)
177 {
178 for (label i = 1; i < dynVerts.size(); i++)
179 {
180 const edge e(dynVerts[i-1], dynVerts[i]);
181 dynEdges.append(e);
182
183 dynUsedPoints[e[0]] = e[0];
184 dynUsedPoints[e[1]] = e[1];
185 }
186 }
187 }
188 }
189
190 // Cull unused points
191 label nUsed = 0;
192 forAll(dynPoints, pointi)
193 {
194 if (dynUsedPoints[pointi] >= 0)
195 {
196 if (nUsed != pointi)
197 {
198 dynPoints[nUsed] = std::move(dynPoints[pointi]);
199 dynUsedPoints[pointi] = nUsed; // The new list location
200 }
201 ++nUsed;
202 }
203 }
204
205 dynPoints.setSize(nUsed);
206
207 // Transfer to normal lists
208 storedPoints().transfer(dynPoints);
209
210 // Renumber edge vertices
211 if (nUsed != dynUsedPoints.size())
212 {
213 for (edge& e : dynEdges)
214 {
215 e[0] = dynUsedPoints[e[0]];
216 e[1] = dynUsedPoints[e[1]];
217 }
218 }
219 storedEdges().transfer(dynEdges);
220
221 return true;
222}
223
224
226(
227 const fileName& filename,
228 const edgeMesh& mesh,
229 IOstreamOption streamOpt,
230 const dictionary&
231)
232{
233 // ASCII only, allow output compression
234 streamOpt.format(IOstream::ASCII);
235
236 const pointField& pointLst = mesh.points();
237 const edgeList& edgeLst = mesh.edges();
238
239 OFstream os(filename, streamOpt);
240 if (!os.good())
241 {
243 << "Cannot write file " << filename << nl
244 << exit(FatalError);
245 }
246
247
248 os << "# Wavefront OBJ file written " << clock::dateTime().c_str() << nl
249 << "o " << os.name().nameLessExt() << nl
250 << nl
251 << "# points : " << pointLst.size() << nl
252 << "# lines : " << edgeLst.size() << nl;
253
254 os << nl
255 << "# <points count=\"" << pointLst.size() << "\">" << nl;
256
257 // Write vertex coords
258 for (const point& p : pointLst)
259 {
260 os << "v " << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
261 }
262
263 os << "# </points>" << nl
264 << nl
265 << "# <edges count=\"" << edgeLst.size() << "\">" << endl;
266
267 // Write line connectivity
268 for (const edge& e : edgeLst)
269 {
270 os << "l " << (e[0] + 1) << " " << (e[1] + 1) << nl;
271 }
272 os << "# </edges>" << endl;
273}
274
275
276// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
void transfer(List< T > &list)
Transfer contents of the argument List into this.
Definition: DynamicListI.H:467
void setSize(const label n)
Same as resize()
Definition: DynamicList.H:221
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
Input from file stream, using an ISstream.
Definition: IFstream.H:57
The IOstreamOption is a simple container for options an IOstream can normally have.
streamFormat format() const noexcept
Get the current stream format.
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
Output to file stream, using an OSstream.
Definition: OFstream.H:57
virtual bool read()
Re-read model coefficients if they have changed.
Sub-ranges of a string with a structure similar to std::match_results, but without the underlying reg...
Definition: SubStrings.H:54
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static std::string dateTime()
Definition: clock.C:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Area discretisation.
Definition: edgeFaMesh.H:56
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Provide a means of reading/writing Alias/Wavefront OBJ format.
Definition: OBJedgeFormat.H:60
virtual bool read(const fileName &name)
Read from file.
Definition: OBJedgeFormat.C:96
A class for handling file names.
Definition: fileName.H:76
virtual bool write()
Write the output fields.
virtual void validate()
Validate the turbulence fields after construction.
Definition: kkLOmega.C:597
A line primitive.
Definition: line.H:68
static constexpr uint64_t npos
Out of range position or size.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
patchWriters clear()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
Foam::SubStrings< StringType > splitSpace(const StringType &str)
Split string into sub-strings at whitespace (TAB, NL, VT, FF, CR, SPC)
Namespace for OpenFOAM.
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition: label.H:66
static label readObjVertices(const SubStrings< string > &tokens, DynamicList< label > &verts)
Definition: OBJedgeFormat.C:46
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333