Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / edgeMesh / edgeFormats / obj / OBJedgeFormat.C
bloba349e71c5547f380e7c95457e177df5e18d4278d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2009-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by
13     the Free Software Foundation, either version 3 of the License, or
14     (at your option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "OBJedgeFormat.H"
27 #include "clock.H"
28 #include "IFstream.H"
29 #include "IStringStream.H"
30 #include "Ostream.H"
31 #include "OFstream.H"
32 #include "ListOps.H"
34 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
36 Foam::fileFormats::OBJedgeFormat::OBJedgeFormat
38     const fileName& filename
41     read(filename);
45 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
47 bool Foam::fileFormats::OBJedgeFormat::read
49     const fileName& filename
52     clear();
54     IFstream is(filename);
55     if (!is.good())
56     {
57         FatalErrorIn
58         (
59             "fileFormats::OBJedgeFormat::read(const fileName&)"
60         )
61             << "Cannot read file " << filename
62             << exit(FatalError);
63     }
65     DynamicList<point> dynPoints;
66     DynamicList<edge> dynEdges;
67     DynamicList<label> dynUsedPoints;
69     while (is.good())
70     {
71         string line = this->getLineNoComment(is);
73         // handle continuations
74         if (line[line.size()-1] == '\\')
75         {
76             line.substr(0, line.size()-1);
77             line += this->getLineNoComment(is);
78         }
80         // Read first word
81         IStringStream lineStream(line);
82         word cmd;
83         lineStream >> cmd;
85         if (cmd == "v")
86         {
87             scalar x, y, z;
88             lineStream >> x >> y >> z;
90             dynPoints.append(point(x, y, z));
91             dynUsedPoints.append(-1);
92         }
93         else if (cmd == "l")
94         {
95             edge edgeRead;
97             // Assume 'l' is followed by space.
98             string::size_type endNum = 1;
100             int nVerts = 0;
101             for (int count = 0; count < 2; ++count)
102             {
103                 string::size_type startNum =
104                     line.find_first_not_of(' ', endNum);
106                 if (startNum == string::npos)
107                 {
108                     break;
109                 }
111                 endNum = line.find(' ', startNum);
113                 string vertexSpec;
114                 if (endNum != string::npos)
115                 {
116                     vertexSpec = line.substr(startNum, endNum-startNum);
117                 }
118                 else
119                 {
120                     vertexSpec = line.substr(startNum, line.size() - startNum);
121                 }
123                 string::size_type slashPos = vertexSpec.find('/');
125                 label vertI = 0;
126                 if (slashPos != string::npos)
127                 {
128                     IStringStream intStream(vertexSpec.substr(0, slashPos));
130                     intStream >> vertI;
131                 }
132                 else
133                 {
134                     IStringStream intStream(vertexSpec);
136                     intStream >> vertI;
137                 }
139                 edgeRead[nVerts++] = (vertI - 1); // change to zero-offset
140             }
142             if (nVerts >= 2)
143             {
144                 dynUsedPoints[edgeRead[0]] = edgeRead[0];
145                 dynUsedPoints[edgeRead[1]] = edgeRead[1];
147                 dynEdges.append(edgeRead);
148             }
149         }
150     }
152     // cull unused points
153     label nUsed = 0;
155     forAll(dynPoints, pointI)
156     {
157         if (dynUsedPoints[pointI] >= 0)
158         {
159             if (nUsed != pointI)
160             {
161                 dynPoints[nUsed] = dynPoints[pointI];
162                 dynUsedPoints[pointI] = nUsed;   // new position
163             }
164             ++nUsed;
165         }
166     }
168     dynPoints.setSize(nUsed);
170     // transfer to normal lists
171     storedPoints().transfer(dynPoints);
173     // renumber edge vertices
174     if (nUsed != dynUsedPoints.size())
175     {
176         forAll(dynEdges, edgeI)
177         {
178             edge& e = dynEdges[edgeI];
180             e[0] = dynUsedPoints[e[0]];
181             e[1] = dynUsedPoints[e[1]];
182         }
183     }
184     storedEdges().transfer(dynEdges);
186     return true;
190 void Foam::fileFormats::OBJedgeFormat::write
192     const fileName& filename,
193     const edgeMesh& mesh
196     const pointField& pointLst = mesh.points();
197     const edgeList& edgeLst = mesh.edges();
199     OFstream os(filename);
200     if (!os.good())
201     {
202         FatalErrorIn
203         (
204             "fileFormats::OBJedgeFormat::write"
205             "(const fileName&, const edgeMesh&)"
206         )
207             << "Cannot open file for writing " << filename
208             << exit(FatalError);
209     }
212     os  << "# Wavefront OBJ file written " << clock::dateTime().c_str() << nl
213         << "o " << os.name().lessExt().name() << nl
214         << nl
215         << "# points : " << pointLst.size() << nl
216         << "# lines  : " << edgeLst.size() << nl;
218     os  << nl
219         << "# <points count=\"" << pointLst.size() << "\">" << nl;
221     // Write vertex coords
222     forAll(pointLst, ptI)
223     {
224         const point& p = pointLst[ptI];
226         os  << "v " << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
227     }
229     os  << "# </points>" << nl
230         << nl
231         << "# <edges count=\"" << edgeLst.size() << "\">" << endl;
233     // Write line connectivity
234     forAll(edgeLst, edgeI)
235     {
236         const edge& e = edgeLst[edgeI];
238         os << "l " << (e[0] + 1) << " " << (e[1] + 1) << nl;
239     }
240     os << "# </edges>" << endl;
244 // ************************************************************************* //