Merge remote-tracking branch 'origin/nr/multiSolverFix' into nextRelease
[foam-extend-3.2.git] / src / surfMesh / surfaceFormats / nas / NASsurfaceFormat.C
blob21669054d416abfbce68bc303b5e8137f0db833b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 "NASsurfaceFormat.H"
28 #include "IFstream.H"
29 #include "IStringStream.H"
31 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
33 template<class Face>
34 Foam::fileFormats::NASsurfaceFormat<Face>::NASsurfaceFormat
36     const fileName& filename
39     this->read(filename);
43 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
45 template<class Face>
46 bool Foam::fileFormats::NASsurfaceFormat<Face>::read
48     const fileName& filename
51     const bool mustTriangulate = this->isTri();
52     this->clear();
54     IFstream is(filename);
55     if (!is.good())
56     {
57         FatalErrorIn
58         (
59             "fileFormats::NASsurfaceFormat::read(const fileName&)"
60         )
61             << "Cannot read file " << filename
62             << exit(FatalError);
63     }
65     // Nastran index of points
66     DynamicList<label>  pointId;
67     DynamicList<point>  dynPoints;
68     DynamicList<Face>   dynFaces;
69     DynamicList<label>  dynZones;
70     DynamicList<label>  dynSizes;
71     Map<label>          lookup;
73     // assume the types are not intermixed
74     // leave faces that didn't have a group in 0
75     bool sorted = true;
76     label zoneI = 0;
78     // Name for face group
79     Map<word> nameLookup;
81     // Ansa tags. Denoted by $ANSA_NAME.
82     // These will appear just before the first use of a type.
83     // We read them and store the PSHELL types which are used to name
84     // the zones.
85     label ansaId = -1;
86     word  ansaType, ansaName;
88     // A single warning per unrecognized command
89     HashSet<word> unhandledCmd;
91     while (is.good())
92     {
93         string line;
94         is.getLine(line);
96         // Ansa extension
97         if (line.substr(0, 10) == "$ANSA_NAME")
98         {
99             string::size_type sem0 = line.find (';', 0);
100             string::size_type sem1 = line.find (';', sem0+1);
101             string::size_type sem2 = line.find (';', sem1+1);
103             if
104             (
105                 sem0 != string::npos
106              && sem1 != string::npos
107              && sem2 != string::npos
108             )
109             {
110                 ansaId = readLabel
111                 (
112                     IStringStream(line.substr(sem0+1, sem1-sem0-1))()
113                 );
114                 ansaType = line.substr(sem1+1, sem2-sem1-1);
116                 string rawName;
117                 is.getLine(rawName);
118                 if (rawName[rawName.size()-1] == '\r')
119                 {
120                     rawName = rawName.substr(1, rawName.size()-2);
121                 }
122                 else
123                 {
124                     rawName = rawName.substr(1, rawName.size()-1);
125                 }
127                 string::stripInvalid<word>(rawName);
128                 ansaName = rawName;
130                 // Info<< "ANSA tag for NastranID:" << ansaId
131                 //     << " of type " << ansaType
132                 //     << " name " << ansaName << endl;
133             }
134         }
137         // Hypermesh extension
138         // $HMNAME COMP                   1"partName"
139         if
140         (
141             line.substr(0, 12) == "$HMNAME COMP"
142          && line.find ('"') != string::npos
143         )
144         {
145             label groupId = readLabel
146             (
147                 IStringStream(line.substr(16, 16))()
148             );
150             IStringStream lineStream(line.substr(32));
152             string rawName;
153             lineStream >> rawName;
154             string::stripInvalid<word>(rawName);
156             word groupName(rawName);
157             nameLookup.insert(groupId, groupName);
159             // Info<< "group " << groupId << " => " << groupName << endl;
160         }
163         // Skip empty or comment
164         if (line.empty() || line[0] == '$')
165         {
166             continue;
167         }
169         // Check if character 72 is continuation
170         if (line.size() > 72 && line[72] == '+')
171         {
172             line = line.substr(0, 72);
174             while (true)
175             {
176                 string buf;
177                 is.getLine(buf);
179                 if (buf.size() > 72 && buf[72] == '+')
180                 {
181                     line += buf.substr(8, 64);
182                 }
183                 else
184                 {
185                     line += buf.substr(8, buf.size()-8);
186                     break;
187                 }
188             }
189         }
192         // Read first word
193         IStringStream lineStream(line);
194         word cmd;
195         lineStream >> cmd;
197         if (cmd == "CTRIA3")
198         {
199             triFace fTri;
201             label groupId = readLabel(IStringStream(line.substr(16,8))());
202             fTri[0] = readLabel(IStringStream(line.substr(24,8))());
203             fTri[1] = readLabel(IStringStream(line.substr(32,8))());
204             fTri[2] = readLabel(IStringStream(line.substr(40,8))());
206             // Convert groupID into zoneId
207             Map<label>::const_iterator fnd = lookup.find(groupId);
208             if (fnd != lookup.end())
209             {
210                 if (zoneI != fnd())
211                 {
212                     // pshell types are intermixed
213                     sorted = false;
214                 }
215                 zoneI = fnd();
216             }
217             else
218             {
219                 zoneI = dynSizes.size();
220                 lookup.insert(groupId, zoneI);
221                 dynSizes.append(0);
222                 // Info<< "zone" << zoneI << " => group " << groupId <<endl;
223             }
225             dynFaces.append(fTri);
226             dynZones.append(zoneI);
227             dynSizes[zoneI]++;
228         }
229         else if (cmd == "CQUAD4")
230         {
231             face fQuad(4);
232             UList<label>& f = static_cast<UList<label>&>(fQuad);
234             label groupId = readLabel(IStringStream(line.substr(16,8))());
235             fQuad[0] = readLabel(IStringStream(line.substr(24,8))());
236             fQuad[1] = readLabel(IStringStream(line.substr(32,8))());
237             fQuad[2] = readLabel(IStringStream(line.substr(40,8))());
238             fQuad[3] = readLabel(IStringStream(line.substr(48,8))());
240             // Convert groupID into zoneId
241             Map<label>::const_iterator fnd = lookup.find(groupId);
242             if (fnd != lookup.end())
243             {
244                 if (zoneI != fnd())
245                 {
246                     // pshell types are intermixed
247                     sorted = false;
248                 }
249                 zoneI = fnd();
250             }
251             else
252             {
253                 zoneI = dynSizes.size();
254                 lookup.insert(groupId, zoneI);
255                 dynSizes.append(0);
256                 // Info<< "zone" << zoneI << " => group " << groupId <<endl;
257             }
260             if (mustTriangulate)
261             {
262                 dynFaces.append(triFace(f[0], f[1], f[2]));
263                 dynFaces.append(triFace(f[0], f[2], f[3]));
264                 dynZones.append(zoneI);
265                 dynZones.append(zoneI);
266                 dynSizes[zoneI] += 2;
267             }
268             else
269             {
270                 dynFaces.append(Face(f));
271                 dynZones.append(zoneI);
272                 dynSizes[zoneI]++;
273             }
274         }
275         else if (cmd == "GRID")
276         {
277             label index = readLabel(IStringStream(line.substr(8,8))());
278             scalar x = parseNASCoord(line.substr(24, 8));
279             scalar y = parseNASCoord(line.substr(32, 8));
280             scalar z = parseNASCoord(line.substr(40, 8));
282             pointId.append(index);
283             dynPoints.append(point(x, y, z));
284         }
285         else if (cmd == "GRID*")
286         {
287             // Long format is on two lines with '*' continuation symbol
288             // on start of second line.
289             // Typical line (spaces compacted)
290             // GRID*      126   0 -5.55999875E+02 -5.68730474E+02
291             // *         2.14897901E+02
293             label index = readLabel(IStringStream(line.substr(8,16))());
294             scalar x = parseNASCoord(line.substr(40, 16));
295             scalar y = parseNASCoord(line.substr(56, 16));
297             is.getLine(line);
298             if (line[0] != '*')
299             {
300                 FatalErrorIn
301                 (
302                     "fileFormats::NASsurfaceFormat::read(const fileName&)"
303                 )
304                     << "Expected continuation symbol '*' when reading GRID*"
305                     << " (double precision coordinate) format" << nl
306                     << "Read:" << line << nl
307                     << "File:" << is.name() << " line:" << is.lineNumber()
308                     << exit(FatalError);
309             }
310             scalar z = parseNASCoord(line.substr(8, 16));
312             pointId.append(index);
313             dynPoints.append(point(x, y, z));
314         }
315         else if (cmd == "PSHELL")
316         {
317             // pshell type for zone names with the Ansa extension
318             label groupId = readLabel(IStringStream(line.substr(8,8))());
320             if (groupId == ansaId && ansaType == "PSHELL")
321             {
322                 nameLookup.insert(ansaId, ansaName);
323                 // Info<< "group " << groupId << " => " << ansaName << endl;
324             }
325         }
326         else if (unhandledCmd.insert(cmd))
327         {
328             Info<< "Unhandled Nastran command " << line << nl
329                 << "File:" << is.name() << " line:" << is.lineNumber()
330                 << endl;
331         }
332     }
334     //    Info<< "Read faces:" << dynFaces.size()
335     //        << " points:" << dynPoints.size()
336     //        << endl;
338     // transfer to normal lists
339     this->storedPoints().transfer(dynPoints);
341     pointId.shrink();
342     dynFaces.shrink();
344     // Build inverse mapping (NASTRAN pointId -> index)
345     Map<label> mapPointId(2*pointId.size());
346     forAll(pointId, i)
347     {
348         mapPointId.insert(pointId[i], i);
349     }
351     // Relabel faces
352     // ~~~~~~~~~~~~~
353     forAll(dynFaces, i)
354     {
355         Face& f = dynFaces[i];
356         forAll(f, fp)
357         {
358             f[fp] = mapPointId[f[fp]];
359         }
360     }
361     pointId.clearStorage();
362     mapPointId.clear();
365     // create default zone names, or from ANSA/Hypermesh information
366     List<word> names(dynSizes.size());
367     forAllConstIter(Map<label>, lookup, iter)
368     {
369         const label zoneI  = iter();
370         const label groupI = iter.key();
372         Map<word>::const_iterator fnd = nameLookup.find(groupI);
373         if (fnd != nameLookup.end())
374         {
375             names[zoneI] = fnd();
376         }
377         else
378         {
379             names[zoneI] = word("zone") + ::Foam::name(zoneI);
380         }
381     }
383     this->sortFacesAndStore(dynFaces.xfer(), dynZones.xfer(), sorted);
385     // add zones, culling empty ones
386     this->addZones(dynSizes, names, true);
388     return true;
392 // ************************************************************************* //