Loading...
Searching...
No Matches
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-2025 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
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& 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 std::string vrtSpec(tok.str());
64
65 if
66 (
67 const auto slash = vrtSpec.find('/');
68 slash != std::string::npos
69 )
70 {
71 vrtSpec.erase(slash);
72 }
73
74 label vertId = readLabel(vrtSpec);
75
76 verts.push_back(vertId - 1);
77 }
78
79 return verts.size();
80}
82} // End namespace Foam
83
84
85
86// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87
89(
90 const fileName& filename
91)
93 read(filename);
94}
95
96
97// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98
100{
101 clear();
102
103 IFstream is(filename);
104 if (!is.good())
105 {
107 << "Cannot read file " << filename
108 << exit(FatalError);
109 }
110
111 DynamicList<point> dynPoints;
112 DynamicList<label> dynVerts;
113 DynamicList<edge> dynEdges;
114 DynamicList<label> dynUsedPoints;
115
116 while (is.good())
117 {
118 string line = this->getLineNoComment(is);
119
120 // Line continuations
121 while (line.removeEnd('\\'))
122 {
123 line += this->getLineNoComment(is);
124 }
125
126 const auto tokens = stringOps::splitSpace(line);
127
128 // Require command and some arguments
129 if (tokens.size() < 2)
130 {
131 continue;
132 }
133
134 const word cmd = word::validate(tokens[0]);
135
136 if (cmd == "v")
137 {
138 // Vertex
139 // v x y z
140
141 dynPoints.append
142 (
143 point
144 (
145 readScalar(tokens[1]),
146 readScalar(tokens[2]),
147 readScalar(tokens[3])
148 )
149 );
150
151 dynUsedPoints.append(-1);
152 }
153 else if (cmd == "l")
154 {
155 // Line
156 // l v1 v2 v3 ...
157 // OR
158 // l v1/vt1 v2/vt2 v3/vt3 ...
159
160 readObjVertices(tokens, dynVerts);
161
162 for (label i = 1; i < dynVerts.size(); i++)
163 {
164 const edge e(dynVerts[i-1], dynVerts[i]);
165 dynEdges.append(e);
166
167 dynUsedPoints[e[0]] = e[0];
168 dynUsedPoints[e[1]] = e[1];
169 }
170 }
171 else if (cmd == "f")
172 {
173 // Face - support for faces with 2 vertices
174
175 // f v1 v2 v3 ...
176 // OR
177 // f v1/vt1 v2/vt2 v3/vt3 ...
178
179 if (readObjVertices(tokens, dynVerts) == 2)
180 {
181 for (label i = 1; i < dynVerts.size(); i++)
182 {
183 const edge e(dynVerts[i-1], dynVerts[i]);
184 dynEdges.append(e);
185
186 dynUsedPoints[e[0]] = e[0];
187 dynUsedPoints[e[1]] = e[1];
188 }
189 }
190 }
191 }
192
193 // Cull unused points
194 label nUsed = 0;
195 forAll(dynPoints, pointi)
196 {
197 if (dynUsedPoints[pointi] >= 0)
198 {
199 if (nUsed != pointi)
200 {
201 dynPoints[nUsed] = std::move(dynPoints[pointi]);
202 dynUsedPoints[pointi] = nUsed; // The new list location
203 }
204 ++nUsed;
205 }
206 }
207
208 dynPoints.setSize(nUsed);
209
210 // Transfer to normal lists
211 storedPoints().transfer(dynPoints);
212
213 // Renumber edge vertices
214 if (nUsed != dynUsedPoints.size())
215 {
216 for (edge& e : dynEdges)
217 {
218 e = edge(dynUsedPoints, e);
219 }
221 storedEdges().transfer(dynEdges);
222
223 return true;
224}
225
226
228(
229 const fileName& filename,
230 const edgeMesh& mesh,
231 IOstreamOption streamOpt,
232 const dictionary&
233)
234{
235 // ASCII only, allow output compression
236 streamOpt.format(IOstreamOption::ASCII);
237
238 const pointField& pointLst = mesh.points();
239 const edgeList& edgeLst = mesh.edges();
240
241 OFstream os(filename, streamOpt);
242 if (!os.good())
243 {
245 << "Cannot write file " << filename << nl
246 << exit(FatalError);
247 }
248
249
250 os << "# Wavefront OBJ file written " << clock::dateTime().c_str() << nl
251 << "o " << os.name().stem() << nl
252 << nl
253 << "# points : " << pointLst.size() << nl
254 << "# lines : " << edgeLst.size() << nl;
255
256 os << nl
257 << "# <points count=\"" << pointLst.size() << "\">" << nl;
258
259 // Write vertex coords
260 for (const point& p : pointLst)
261 {
262 os << "v " << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
263 }
264
265 os << "# </points>" << nl
266 << nl
267 << "# <edges count=\"" << edgeLst.size() << "\">" << endl;
268
269 // Write line connectivity
270 for (const edge& e : edgeLst)
271 {
272 os << "l " << (e[0] + 1) << " " << (e[1] + 1) << nl;
273 }
274 os << "# </edges>" << endl;
275}
276
277
278// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void setSize(const label n)
Alias for resize().
void append(const T &val)
Copy append an element to the end of this list.
void push_back(const T &val)
Copy append an element to the end of this list.
Input from file stream as an ISstream, normally using std::ifstream for the actual input.
Definition IFstream.H:55
A simple container for options an IOstream can normally have.
streamFormat format() const noexcept
Get the current stream format.
@ ASCII
"ascii" (normal default)
bool good() const noexcept
True if next operation might succeed.
Definition IOstream.H:281
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
Sub-ranges of a string with a structure similar to std::match_results, but without the underlying reg...
Definition SubStrings.H:49
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static std::string dateTime()
The current wall-clock date/time (in local time) as a string in ISO-8601 format (yyyy-mm-ddThh:mm:ss)...
Definition clock.C:53
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Area discretisation.
Definition edgeFaMesh.H:50
pointField & storedPoints() noexcept
Non-const access to global points.
Definition edgeMeshI.H:24
edgeMesh(const faMesh &mesh)
Construct finite-area edge mesh faMesh reference.
Definition edgeFaMesh.H:58
edgeList & storedEdges() noexcept
Non-const access to the edges.
Definition edgeMeshI.H:30
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
OBJedgeFormat(const fileName &filename)
Construct from file name.
virtual bool read(const fileName &filename) override
Read from file.
static void write(const fileName &filename, const edgeMesh &mesh, IOstreamOption=IOstreamOption(), const dictionary &options=dictionary::null)
Write edge mesh to file in OBJ format.
static string getLineNoComment(ISstream &is, const char comment='#')
Read non-comment line.
A class for handling file names.
Definition fileName.H:75
A line primitive.
Definition line.H:180
A class for handling words, derived from Foam::string.
Definition word.H:66
static word validate(const std::string &s, const bool prefix=false)
Construct validated word (no invalid characters).
Definition word.C:39
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
surface1 clear()
Foam::SubStrings splitSpace(const std::string &str, std::string::size_type pos=0)
Split string into sub-strings at whitespace (TAB, NL, VT, FF, CR, SPC).
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition label.H:63
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
static label readObjVertices(const SubStrings &tokens, DynamicList< label > &verts)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299