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 -------------------------------------------------------------------------------
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 "OBJedgeFormat.H"
30 #include "clock.H"
31 #include "Fstream.H"
32 #include "Ostream.H"
33 #include "stringOps.H"
34 
35 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
36 
37 namespace 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 ...
45 static label readObjVertices
46 (
47  const SubStrings<string>& tokens,
48  DynamicList<label>& verts
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 
85 Foam::fileFormats::OBJedgeFormat::OBJedgeFormat
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 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fileFormats::OBJedgeFormat::read
virtual bool read(const fileName &name)
Read from file.
Definition: OBJedgeFormat.C:96
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
OBJedgeFormat.H
Foam::IFstream
Input from file stream, using an ISstream.
Definition: IFstream.H:53
Foam::DynamicList< label >
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::IOstreamOption::format
streamFormat format() const noexcept
Get the current stream format.
Definition: IOstreamOption.H:286
Foam::word::validate
static word validate(const std::string &s, const bool prefix=false)
Construct validated word (no invalid characters).
Definition: word.C:45
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::readObjVertices
static label readObjVertices(const SubStrings< string > &tokens, DynamicList< label > &verts)
Definition: OBJedgeFormat.C:46
Foam::stringOps::splitSpace
Foam::SubStrings< StringType > splitSpace(const StringType &str)
Split string into sub-strings at whitespace (TAB, NL, VT, FF, CR, SPC)
Definition: stringOpsTemplates.C:246
Foam::SubStrings
Sub-ranges of a string with a structure similar to std::match_results, but without the underlying reg...
Definition: CStringList.H:63
Foam::DynamicList::clear
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
Foam::edgeMesh::storedPoints
pointField & storedPoints() noexcept
Non-const access to global points.
Definition: edgeMeshI.H:31
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::IOstream::good
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::Field< vector >
Foam::fileFormats::OBJedgeFormat::write
static void write(const fileName &filename, const edgeMesh &mesh, IOstreamOption=IOstreamOption(), const dictionary &options=dictionary::null)
Write edge mesh to file.
Definition: OBJedgeFormat.C:226
Foam::clock::dateTime
static std::string dateTime()
Definition: clock.C:60
clock.H
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
Foam::DynamicList::setSize
void setSize(const label n)
Same as resize()
Definition: DynamicList.H:224
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::edgeMesh::clear
virtual void clear()
Clear all storage.
Definition: edgeMesh.C:116
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Ostream.H
Foam::IOstreamOption::ASCII
"ascii" (normal default)
Definition: IOstreamOption.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Fstream.H
Foam::Vector< scalar >
Foam::readLabel
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition: label.H:66
Foam::List< edge >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::line
A line primitive.
Definition: line.H:53
Foam::edgeMesh::storedEdges
edgeList & storedEdges() noexcept
Non-const access to the edges.
Definition: edgeMeshI.H:37
stringOps.H
Foam::edgeMesh
Mesh data needed to do the Finite Area discretisation.
Definition: edgeFaMesh.H:53