initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / surfMesh / surfaceFormats / tri / TRIsurfaceFormatCore.C
bloba6b9c999120b8119f52c06a73879cedfae4a5d08
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "TRIsurfaceFormat.H"
28 #include "IFstream.H"
29 #include "IOmanip.H"
30 #include "IStringStream.H"
31 #include "Map.H"
33 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
35 Foam::fileFormats::TRIsurfaceFormatCore::TRIsurfaceFormatCore
37     const fileName& filename
40     sorted_(true),
41     points_(0),
42     zoneIds_(0),
43     sizes_(0)
45     read(filename);
49 // * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
51 Foam::fileFormats::TRIsurfaceFormatCore::~TRIsurfaceFormatCore()
55 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
57 bool Foam::fileFormats::TRIsurfaceFormatCore::read
59     const fileName& filename
62     this->clear();
63     sorted_ = true;
65     IFstream is(filename);
66     if (!is.good())
67     {
68         FatalErrorIn
69         (
70             "fileFormats::TRIsurfaceFormatCore::read(const fileName&)"
71         )
72             << "Cannot read file " << filename
73             << exit(FatalError);
74     }
76     // uses similar structure as STL, just some points
77     // the rest of the reader resembles the STL binary reader
78     DynamicList<point> dynPoints;
79     DynamicList<label> dynZones;
80     DynamicList<label> dynSizes;
81     HashTable<label>   lookup;
83     // place faces without a group in zone0
84     label zoneI = 0;
85     dynSizes.append(zoneI);
86     lookup.insert("zoneI", zoneI);
88     while (is.good())
89     {
90         string line = this->getLineNoComment(is);
92         // handle continuations ?
93         //          if (line[line.size()-1] == '\\')
94         //          {
95         //              line.substr(0, line.size()-1);
96         //              line += this->getLineNoComment(is);
97         //          }
99         IStringStream lineStream(line);
101         point p
102         (
103             readScalar(lineStream),
104             readScalar(lineStream),
105             readScalar(lineStream)
106         );
108         if (!lineStream) break;
110         dynPoints.append(p);
111         dynPoints.append
112         (
113             point
114             (
115                 readScalar(lineStream),
116                 readScalar(lineStream),
117                 readScalar(lineStream)
118             )
119         );
120         dynPoints.append
121         (
122             point
123             (
124                 readScalar(lineStream),
125                 readScalar(lineStream),
126                 readScalar(lineStream)
127             )
128         );
130         // zone/colour in .tri file starts with 0x. Skip.
131         // ie, instead of having 0xFF, skip 0 and leave xFF to
132         // get read as a word and name it "zoneFF"
134         char zero;
135         lineStream >> zero;
137         word rawName(lineStream);
138         word name("zone" + rawName(1, rawName.size()-1));
140         HashTable<label>::const_iterator fnd = lookup.find(name);
141         if (fnd != lookup.end())
142         {
143             if (zoneI != fnd())
144             {
145                 // group appeared out of order
146                 sorted_ = false;
147             }
148             zoneI = fnd();
149         }
150         else
151         {
152             zoneI = dynSizes.size();
153             lookup.insert(name, zoneI);
154             dynSizes.append(0);
155         }
157         dynZones.append(zoneI);
158         dynSizes[zoneI]++;
159     }
161     // skip empty groups
162     label nZone = 0;
163     forAll(dynSizes, zoneI)
164     {
165         if (dynSizes[zoneI])
166         {
167             if (nZone != zoneI)
168             {
169                 dynSizes[nZone] = dynSizes[zoneI];
170             }
171             nZone++;
172         }
173     }
174     // truncate addressed size
175     dynSizes.setCapacity(nZone);
177     // transfer to normal lists
178     points_.transfer(dynPoints);
179     zoneIds_.transfer(dynZones);
180     sizes_.transfer(dynSizes);
182     return true;
186 // ************************************************************************* //