File: objToVTK.C

package info (click to toggle)
openfoam 4.1%2Bdfsg1-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 163,028 kB
  • ctags: 58,990
  • sloc: cpp: 830,760; sh: 10,227; ansic: 8,215; xml: 745; lex: 437; awk: 194; sed: 91; makefile: 77; python: 18
file content (304 lines) | stat: -rw-r--r-- 6,954 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-2015 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

Application
    objToVTK

Description
    Read obj line (not surface!) file and convert into vtk.

\*---------------------------------------------------------------------------*/

#include "argList.H"
#include "OFstream.H"
#include <fstream>
#include <sstream>
#include "IStringStream.H"
#include "point.H"
#include "DynamicList.H"


using namespace Foam;


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

string getLine(std::ifstream& is)
{
    string line;
    do
    {
        std::getline(is, line);
    }
    while (line.size() && line[0] == '#');

    return line;
}


// Read space-separated vertices (with optional '/' arguments)
labelList parseVertices(const string& line)
{
    DynamicList<label> verts;

    // Assume 'l' is followed by space.
    string::size_type endNum = 1;

    do
    {
        string::size_type startNum = line.find_first_not_of(' ', endNum);

        if (startNum == string::npos)
        {
            break;
        }

        endNum = line.find(' ', startNum);

        string vertexSpec;
        if (endNum != string::npos)
        {
            vertexSpec = line.substr(startNum, endNum-startNum);
        }
        else
        {
            vertexSpec = line.substr(startNum, line.size() - startNum);
        }

        string::size_type slashPos = vertexSpec.find('/');

        label vertI = 0;
        if (slashPos != string::npos)
        {
            IStringStream intStream(vertexSpec.substr(0, slashPos));

            intStream >> vertI;
        }
        else
        {
            IStringStream intStream(vertexSpec);

            intStream >> vertI;
        }
        verts.append(vertI - 1);
    }
    while (true);

    return verts.shrink();
}



int main(int argc, char *argv[])
{
    argList::noParallel();
    argList::validArgs.append("OBJ file");
    argList::validArgs.append("output VTK file");
    argList args(argc, argv);

    const fileName objName = args[1];
    const fileName outName = args[2];

    std::ifstream OBJfile(objName.c_str());

    if (!OBJfile.good())
    {
        FatalErrorInFunction
            << "Cannot read file " << objName << exit(FatalError);
    }

    // Points and lines
    DynamicList<point> points;
    DynamicList<vector> pointNormals;
    DynamicList<labelList> polyLines;
    DynamicList<labelList> polygons;

    bool hasWarned = false;

    label lineNo = 0;
    while (OBJfile.good())
    {
        string line = getLine(OBJfile);
        lineNo++;

        // Read first word
        IStringStream lineStream(line);
        word cmd;
        lineStream >> cmd;

        if (cmd == "v")
        {
            scalar x, y, z;

            lineStream >> x >> y >> z;

            points.append(point(x, y, z));
        }
        else if (cmd == "vn")
        {
            scalar x, y, z;

            lineStream >> x >> y >> z;

            pointNormals.append(vector(x, y, z));
        }
        else if (cmd == "l")
        {
            polyLines.append(parseVertices(line));
        }
        else if (cmd == "f")
        {
            polygons.append(parseVertices(line));
        }
        else if (cmd != "")
        {
            if (!hasWarned)
            {
                hasWarned = true;

                WarningInFunction
                    << "Unrecognized OBJ command " << cmd << nl
                    << "In line " << lineStream.str()
                    << " at linenumber " << lineNo << nl
                    << "Only recognized commands are 'v' and 'l'.\n"
                    << "If this is a surface command use surfaceConvert instead"
                    << " to convert to a file format that can be read by VTK"
                    << endl;
            }
        }
    }


    //
    // Write as vtk 'polydata' file
    //


    OFstream outFile(outName);

    outFile
        << "# vtk DataFile Version 2.0\n"
        << objName << nl
        << "ASCII\n"
        << "DATASET POLYDATA\n"
        << "POINTS " << points.size() << " float\n";

    forAll(points, i)
    {
        const point& pt = points[i];

        outFile << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
    }

    outFile
        << "VERTICES " << points.size() << ' ' << (2 * points.size()) << nl;

    forAll(points, i)
    {
        outFile << 1 << ' ' << i << nl;
    }

    label nItems = 0;
    forAll(polyLines, polyI)
    {
        nItems += polyLines[polyI].size() + 1;
    }

    outFile
        << "LINES " << polyLines.size() << ' ' << nItems << nl;

    forAll(polyLines, polyI)
    {
        const labelList& line = polyLines[polyI];

        outFile << line.size();

        forAll(line, i)
        {
            outFile << ' ' << line[i];
        }
        outFile << nl;
    }


    nItems = 0;
    forAll(polygons, polyI)
    {
        nItems += polygons[polyI].size() + 1;
    }

    outFile
        << "POLYGONS " << polygons.size() << ' ' << nItems << nl;

    forAll(polygons, polyI)
    {
        const labelList& line = polygons[polyI];

        outFile << line.size();

        forAll(line, i)
        {
            outFile << ' ' << line[i];
        }
        outFile << nl;
    }


    outFile
        << "POINT_DATA " << points.size() << nl
        << "SCALARS pointID float 1\n"
        << "LOOKUP_TABLE default\n";

    forAll(points, i)
    {
        outFile << i;

        if ((i % 10) == 1)
        {
            outFile << nl;
        }
        else
        {
            outFile << ' ';
        }
    }

    if (!pointNormals.empty())
    {
        outFile << nl << "NORMALS pointNormals float\n";

        forAll(pointNormals, i)
        {
            const vector& n = pointNormals[i];

            outFile << n.x() << ' ' << n.y() << ' ' << n.z() << nl;
        }
    }

    Info<< "End\n" << endl;

    return 0;
}


// ************************************************************************* //