Updated .L files to assign properly the function 'yywrap()' for Flex 2.6.0 and newer.
[foam-extend-3.2.git] / src / surfMesh / surfaceFormats / nas / NASsurfaceFormat.C
blob480c6f76b15c29077eaf560c6ee5463a844d4221
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "NASsurfaceFormat.H"
27 #include "IFstream.H"
28 #include "IStringStream.H"
30 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
32 template<class Face>
33 Foam::fileFormats::NASsurfaceFormat<Face>::NASsurfaceFormat
35     const fileName& filename
38     this->read(filename);
42 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
44 template<class Face>
45 bool Foam::fileFormats::NASsurfaceFormat<Face>::read
47     const fileName& filename
50     const bool mustTriangulate = this->isTri();
51     this->clear();
53     IFstream is(filename);
54     if (!is.good())
55     {
56         FatalErrorIn
57         (
58             "fileFormats::NASsurfaceFormat::read(const fileName&)"
59         )
60             << "Cannot read file " << filename
61             << exit(FatalError);
62     }
64     // Nastran index of points
65     dynamicLabelList  pointId;
66     DynamicList<point>  dynPoints;
67     DynamicList<Face>   dynFaces;
68     dynamicLabelList  dynZones;
69     dynamicLabelList  dynSizes;
70     Map<label>          lookup;
72     // assume the types are not intermixed
73     // leave faces that didn't have a group in 0
74     bool sorted = true;
75     label zoneI = 0;
77     // Name for face group
78     Map<word> nameLookup;
80     // Ansa tags. Denoted by $ANSA_NAME.
81     // These will appear just before the first use of a type.
82     // We read them and store the PSHELL types which are used to name
83     // the zones.
84     label ansaId = -1;
85     word  ansaType, ansaName;
87     // A single warning per unrecognized command
88     HashSet<word> unhandledCmd;
90     while (is.good())
91     {
92         string line;
93         is.getLine(line);
95         // Ansa extension
96         if (line.substr(0, 10) == "$ANSA_NAME")
97         {
98             string::size_type sem0 = line.find (';', 0);
99             string::size_type sem1 = line.find (';', sem0+1);
100             string::size_type sem2 = line.find (';', sem1+1);
102             if
103             (
104                 sem0 != string::npos
105              && sem1 != string::npos
106              && sem2 != string::npos
107             )
108             {
109                 ansaId = readLabel
110                 (
111                     IStringStream(line.substr(sem0+1, sem1-sem0-1))()
112                 );
113                 ansaType = line.substr(sem1+1, sem2-sem1-1);
115                 string rawName;
116                 is.getLine(rawName);
117                 if (rawName[rawName.size()-1] == '\r')
118                 {
119                     rawName = rawName.substr(1, rawName.size()-2);
120                 }
121                 else
122                 {
123                     rawName = rawName.substr(1, rawName.size()-1);
124                 }
126                 string::stripInvalid<word>(rawName);
127                 ansaName = rawName;
129                 // Info<< "ANSA tag for NastranID:" << ansaId
130                 //     << " of type " << ansaType
131                 //     << " name " << ansaName << endl;
132             }
133         }
136         // Hypermesh extension
137         // $HMNAME COMP                   1"partName"
138         if
139         (
140             line.substr(0, 12) == "$HMNAME COMP"
141          && line.find ('"') != string::npos
142         )
143         {
144             label groupId = readLabel
145             (
146                 IStringStream(line.substr(16, 16))()
147             );
149             IStringStream lineStream(line.substr(32));
151             string rawName;
152             lineStream >> rawName;
153             string::stripInvalid<word>(rawName);
155             word groupName(rawName);
156             nameLookup.insert(groupId, groupName);
158             // Info<< "group " << groupId << " => " << groupName << endl;
159         }
162         // Skip empty or comment
163         if (line.empty() || line[0] == '$')
164         {
165             continue;
166         }
168         // Check if character 72 is continuation
169         if (line.size() > 72 && line[72] == '+')
170         {
171             line = line.substr(0, 72);
173             while (true)
174             {
175                 string buf;
176                 is.getLine(buf);
178                 if (buf.size() > 72 && buf[72] == '+')
179                 {
180                     line += buf.substr(8, 64);
181                 }
182                 else
183                 {
184                     line += buf.substr(8, buf.size()-8);
185                     break;
186                 }
187             }
188         }
191         // Read first word
192         IStringStream lineStream(line);
193         word cmd;
194         lineStream >> cmd;
196         if (cmd == "CTRIA3")
197         {
198             triFace fTri;
200             label groupId = readLabel(IStringStream(line.substr(16,8))());
201             fTri[0] = readLabel(IStringStream(line.substr(24,8))());
202             fTri[1] = readLabel(IStringStream(line.substr(32,8))());
203             fTri[2] = readLabel(IStringStream(line.substr(40,8))());
205             // Convert groupID into zoneId
206             Map<label>::const_iterator fnd = lookup.find(groupId);
207             if (fnd != lookup.end())
208             {
209                 if (zoneI != fnd())
210                 {
211                     // pshell types are intermixed
212                     sorted = false;
213                 }
214                 zoneI = fnd();
215             }
216             else
217             {
218                 zoneI = dynSizes.size();
219                 lookup.insert(groupId, zoneI);
220                 dynSizes.append(0);
221                 // Info<< "zone" << zoneI << " => group " << groupId <<endl;
222             }
224             dynFaces.append(fTri);
225             dynZones.append(zoneI);
226             dynSizes[zoneI]++;
227         }
228         else if (cmd == "CQUAD4")
229         {
230             face fQuad(4);
231             UList<label>& f = static_cast<UList<label>&>(fQuad);
233             label groupId = readLabel(IStringStream(line.substr(16,8))());
234             fQuad[0] = readLabel(IStringStream(line.substr(24,8))());
235             fQuad[1] = readLabel(IStringStream(line.substr(32,8))());
236             fQuad[2] = readLabel(IStringStream(line.substr(40,8))());
237             fQuad[3] = readLabel(IStringStream(line.substr(48,8))());
239             // Convert groupID into zoneId
240             Map<label>::const_iterator fnd = lookup.find(groupId);
241             if (fnd != lookup.end())
242             {
243                 if (zoneI != fnd())
244                 {
245                     // pshell types are intermixed
246                     sorted = false;
247                 }
248                 zoneI = fnd();
249             }
250             else
251             {
252                 zoneI = dynSizes.size();
253                 lookup.insert(groupId, zoneI);
254                 dynSizes.append(0);
255                 // Info<< "zone" << zoneI << " => group " << groupId <<endl;
256             }
259             if (mustTriangulate)
260             {
261                 dynFaces.append(triFace(f[0], f[1], f[2]));
262                 dynFaces.append(triFace(f[0], f[2], f[3]));
263                 dynZones.append(zoneI);
264                 dynZones.append(zoneI);
265                 dynSizes[zoneI] += 2;
266             }
267             else
268             {
269                 dynFaces.append(Face(f));
270                 dynZones.append(zoneI);
271                 dynSizes[zoneI]++;
272             }
273         }
274         else if (cmd == "GRID")
275         {
276             label index = readLabel(IStringStream(line.substr(8,8))());
277             scalar x = parseNASCoord(line.substr(24, 8));
278             scalar y = parseNASCoord(line.substr(32, 8));
279             scalar z = parseNASCoord(line.substr(40, 8));
281             pointId.append(index);
282             dynPoints.append(point(x, y, z));
283         }
284         else if (cmd == "GRID*")
285         {
286             // Long format is on two lines with '*' continuation symbol
287             // on start of second line.
288             // Typical line (spaces compacted)
289             // GRID*      126   0 -5.55999875E+02 -5.68730474E+02
290             // *         2.14897901E+02
292             label index = readLabel(IStringStream(line.substr(8,16))());
293             scalar x = parseNASCoord(line.substr(40, 16));
294             scalar y = parseNASCoord(line.substr(56, 16));
296             is.getLine(line);
297             if (line[0] != '*')
298             {
299                 FatalErrorIn
300                 (
301                     "fileFormats::NASsurfaceFormat::read(const fileName&)"
302                 )
303                     << "Expected continuation symbol '*' when reading GRID*"
304                     << " (double precision coordinate) format" << nl
305                     << "Read:" << line << nl
306                     << "File:" << is.name() << " line:" << is.lineNumber()
307                     << exit(FatalError);
308             }
309             scalar z = parseNASCoord(line.substr(8, 16));
311             pointId.append(index);
312             dynPoints.append(point(x, y, z));
313         }
314         else if (cmd == "PSHELL")
315         {
316             // pshell type for zone names with the Ansa extension
317             label groupId = readLabel(IStringStream(line.substr(8,8))());
319             if (groupId == ansaId && ansaType == "PSHELL")
320             {
321                 nameLookup.insert(ansaId, ansaName);
322                 // Info<< "group " << groupId << " => " << ansaName << endl;
323             }
324         }
325         else if (unhandledCmd.insert(cmd))
326         {
327             Info<< "Unhandled Nastran command " << line << nl
328                 << "File:" << is.name() << " line:" << is.lineNumber()
329                 << endl;
330         }
331     }
333     //    Info<< "Read faces:" << dynFaces.size()
334     //        << " points:" << dynPoints.size()
335     //        << endl;
337     // transfer to normal lists
338     this->storedPoints().transfer(dynPoints);
340     pointId.shrink();
341     dynFaces.shrink();
343     // Build inverse mapping (NASTRAN pointId -> index)
344     Map<label> mapPointId(2*pointId.size());
345     forAll(pointId, i)
346     {
347         mapPointId.insert(pointId[i], i);
348     }
350     // Relabel faces
351     // ~~~~~~~~~~~~~
352     forAll(dynFaces, i)
353     {
354         Face& f = dynFaces[i];
355         forAll(f, fp)
356         {
357             f[fp] = mapPointId[f[fp]];
358         }
359     }
360     pointId.clearStorage();
361     mapPointId.clear();
364     // create default zone names, or from ANSA/Hypermesh information
365     List<word> names(dynSizes.size());
366     forAllConstIter(Map<label>, lookup, iter)
367     {
368         const label zoneI  = iter();
369         const label groupI = iter.key();
371         Map<word>::const_iterator fnd = nameLookup.find(groupI);
372         if (fnd != nameLookup.end())
373         {
374             names[zoneI] = fnd();
375         }
376         else
377         {
378             names[zoneI] = word("zone") + ::Foam::name(zoneI);
379         }
380     }
382     this->sortFacesAndStore(dynFaces.xfer(), dynZones.xfer(), sorted);
384     // add zones, culling empty ones
385     this->addZones(dynSizes, names, true);
387     return true;
391 // ************************************************************************* //