readNAS.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Description
25  Nastran surface reader.
26 
27  - Uses the Ansa "$ANSA_NAME" or the Hypermesh "$HMNAME COMP" extensions
28  to obtain patch names.
29  - Handles Nastran short and long formats, but not free format.
30  - Properly handles the Nastran compact floating point notation: \n
31  \verbatim
32  GRID 28 10.20269-.030265-2.358-8
33  \endverbatim
34 
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "triSurface.H"
39 #include "IFstream.H"
40 #include "IStringStream.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
48 
49 // Do weird things to extract number
50 static scalar parseNASCoord(const string& s)
51 {
52  size_t expSign = s.find_last_of("+-");
53 
54  if (expSign != string::npos && expSign > 0 && !isspace(s[expSign-1]))
55  {
56  scalar mantissa = readScalar(IStringStream(s.substr(0, expSign))());
57  scalar exponent = readScalar(IStringStream(s.substr(expSign+1))());
58 
59  if (s[expSign] == '-')
60  {
61  exponent = -exponent;
62  }
63  return mantissa*pow(10, exponent);
64  }
65  else
66  {
67  return readScalar(IStringStream(s)());
68  }
69 }
70 
71 
72 bool triSurface::readNAS(const fileName& fName)
73 {
74  IFstream is(fName);
75 
76  if (!is.good())
77  {
79  << "Cannot read file " << fName
80  << exit(FatalError);
81  }
82 
83  // coordinates of point
85  // Nastran index of point
86  DynamicList<label> indices;
87  // Faces in terms of Nastran point indices
89  // From face group to patch
90  Map<label> groupToPatch;
91  label nPatches = 0;
92  // Name for face group
93  Map<word> groupToName;
94 
95  // Ansa tags. Denoted by $ANSA_NAME. These will appear just before the
96  // first use of a type. We read them and store the pshell types which
97  // are used to name the patches.
98  label ansaId = -1;
99  word ansaType;
100  string ansaName;
101 
102  // A single warning per unrecognized command
103  HashSet<word> unhandledCmd;
104 
105  while (is.good())
106  {
107  string line;
108  is.getLine(line);
109 
110  // Ansa extension
111  if (line.substr(0, 10) == "$ANSA_NAME")
112  {
113  string::size_type sem0 = line.find (';', 0);
114  string::size_type sem1 = line.find (';', sem0+1);
115  string::size_type sem2 = line.find (';', sem1+1);
116 
117  if
118  (
119  sem0 != string::npos
120  && sem1 != string::npos
121  && sem2 != string::npos
122  )
123  {
124  ansaId = readLabel
125  (
126  IStringStream(line.substr(sem0+1, sem1-sem0-1))()
127  );
128  ansaType = line.substr(sem1+1, sem2-sem1-1);
129 
130  string nameString;
131  is.getLine(ansaName);
132  if (ansaName[ansaName.size()-1] == '\r')
133  {
134  ansaName = ansaName.substr(1, ansaName.size()-2);
135  }
136  else
137  {
138  ansaName = ansaName.substr(1, ansaName.size()-1);
139  }
140 
141  // Info<< "ANSA tag for NastranID:" << ansaId
142  // << " of type " << ansaType
143  // << " name " << ansaName << endl;
144  }
145  }
146 
147 
148  // Hypermesh extension
149  // $HMNAME COMP 1"partName"
150  if
151  (
152  line.substr(0, 12) == "$HMNAME COMP"
153  && line.find ('"') != string::npos
154  )
155  {
156  label groupId = readLabel
157  (
158  IStringStream(line.substr(16, 16))()
159  );
160 
161  IStringStream lineStream(line.substr(32));
162 
163  string rawName;
164  lineStream >> rawName;
165 
166  groupToName.insert(groupId, string::validate<word>(rawName));
167  Info<< "group " << groupId << " => " << rawName << endl;
168  }
169 
170 
171  if (line.empty() || line[0] == '$')
172  {
173  // Skip empty or comment
174  continue;
175  }
176 
177  // Check if character 72 is continuation
178  if (line.size() > 72 && line[72] == '+')
179  {
180  line = line.substr(0, 72);
181 
182  while (true)
183  {
184  string buf;
185  is.getLine(buf);
186 
187  if (buf.size() > 72 && buf[72]=='+')
188  {
189  line += buf.substr(8, 64);
190  }
191  else
192  {
193  line += buf.substr(8, buf.size()-8);
194  break;
195  }
196  }
197  }
198 
199  // Read first word
200  IStringStream lineStream(line);
201  word cmd;
202  lineStream >> cmd;
203 
204  if (cmd == "CTRIA3")
205  {
206  label groupId = readLabel(IStringStream(line.substr(16,8))());
207  label a = readLabel(IStringStream(line.substr(24,8))());
208  label b = readLabel(IStringStream(line.substr(32,8))());
209  label c = readLabel(IStringStream(line.substr(40,8))());
210 
211 
212  // Convert group into patch
213  Map<label>::const_iterator iter = groupToPatch.find(groupId);
214 
215  label patchi;
216  if (iter == groupToPatch.end())
217  {
218  patchi = nPatches++;
219  groupToPatch.insert(groupId, patchi);
220  Info<< "patch " << patchi << " => group " << groupId << endl;
221  }
222  else
223  {
224  patchi = iter();
225  }
226 
227  faces.append(labelledTri(a, b, c, patchi));
228  }
229  else if (cmd == "CQUAD4")
230  {
231  label groupId = readLabel(IStringStream(line.substr(16,8))());
232  label a = readLabel(IStringStream(line.substr(24,8))());
233  label b = readLabel(IStringStream(line.substr(32,8))());
234  label c = readLabel(IStringStream(line.substr(40,8))());
235  label d = readLabel(IStringStream(line.substr(48,8))());
236 
237  // Convert group into patch
238  Map<label>::const_iterator iter = groupToPatch.find(groupId);
239 
240  label patchi;
241  if (iter == groupToPatch.end())
242  {
243  patchi = nPatches++;
244  groupToPatch.insert(groupId, patchi);
245  Info<< "patch " << patchi << " => group " << groupId << endl;
246  }
247  else
248  {
249  patchi = iter();
250  }
251 
252  faces.append(labelledTri(a, b, c, patchi));
253  faces.append(labelledTri(c, d, a, patchi));
254  }
255  else if (cmd == "PSHELL")
256  {
257  // Read shell type since group gives patchnames
258  label groupId = readLabel(IStringStream(line.substr(8,8))());
259  if (groupId == ansaId && ansaType == "PSHELL")
260  {
261  groupToName.insert(groupId, string::validate<word>(ansaName));
262  Info<< "group " << groupId << " => " << ansaName << endl;
263  }
264  }
265  else if (cmd == "GRID")
266  {
267  label index = readLabel(IStringStream(line.substr(8,8))());
268  scalar x = parseNASCoord(line.substr(24, 8));
269  scalar y = parseNASCoord(line.substr(32, 8));
270  scalar z = parseNASCoord(line.substr(40, 8));
271 
272  indices.append(index);
273  points.append(point(x, y, z));
274  }
275  else if (cmd == "GRID*")
276  {
277  // Long format is on two lines with '*' continuation symbol
278  // on start of second line.
279  // Typical line (spaces compacted)
280  // GRID* 126 0 -5.55999875E+02 -5.68730474E+02
281  // * 2.14897901E+02
282 
283  label index = readLabel(IStringStream(line.substr(8,16))());
284  scalar x = parseNASCoord(line.substr(40, 16));
285  scalar y = parseNASCoord(line.substr(56, 16));
286 
287  is.getLine(line);
288  if (line[0] != '*')
289  {
291  << "Expected continuation symbol '*' when reading GRID*"
292  << " (double precision coordinate) output" << nl
293  << "Read:" << line << nl
294  << "File:" << is.name()
295  << " line:" << is.lineNumber()
296  << exit(FatalError);
297  }
298  scalar z = parseNASCoord(line.substr(8, 16));
299 
300  indices.append(index);
301  points.append(point(x, y, z));
302  }
303  else if (unhandledCmd.insert(cmd))
304  {
305  Info<< "Unhandled Nastran command " << line << nl
306  << "File:" << is.name() << " line:" << is.lineNumber() << endl;
307  }
308  }
309 
310  points.shrink();
311  indices.shrink();
312  faces.shrink();
313 
314 
315  Info<< "Read triangles:" << faces.size() << " points:" << points.size()
316  << endl;
317 
318  {
319  // Build inverse mapping (index to point)
320  Map<label> indexToPoint(2*indices.size());
321  forAll(indices, i)
322  {
323  indexToPoint.insert(indices[i], i);
324  }
325 
326  // Relabel faces
327  forAll(faces, i)
328  {
329  labelledTri& f = faces[i];
330 
331  f[0] = indexToPoint[f[0]];
332  f[1] = indexToPoint[f[1]];
333  f[2] = indexToPoint[f[2]];
334  }
335  }
336 
337 
338  // Convert groupToPatch to patchList.
340 
341  forAllConstIter(Map<word>, groupToName, iter)
342  {
343  label patchi = groupToPatch[iter.key()];
344 
345  patches[patchi] = geometricSurfacePatch
346  (
347  "empty",
348  iter(),
349  patchi
350  );
351  }
352 
353  Info<< "patches:" << patches << endl;
354 
355  // Transfer DynamicLists to straight ones.
356  pointField allPoints(points.xfer());
357 
358  // Create triSurface
359  *this = triSurface(faces, patches, allPoints, true);
360 
361  return true;
362 }
363 
364 
365 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
366 
367 } // End namespace Foam
368 
369 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:402
A HashTable with keys but without contents.
Definition: HashSet.H:59
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A line primitive.
Definition: line.H:56
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:106
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Field< point > & points() const
Return reference to global points.
const fileName & name() const
Return the name of the stream.
Definition: IFstream.H:116
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
bool insert(const label &, const T &newElmt)
Insert a new hashedEntry.
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:138
scalar y
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from string.
Definition: word.H:59
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:314
Xfer< List< T > > xfer()
Transfer contents to the Xfer container as a plain List.
Definition: DynamicListI.H:283
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
Triangle with additional region number.
Definition: labelledTri.H:57
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:240
static const char nl
Definition: Ostream.H:262
Input from file stream.
Definition: IFstream.H:81
labelList f(nPoints)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
triSurface()
Construct null.
Definition: triSurface.C:594
label patchi
vector point
Point is a vector.
Definition: point.H:41
bool isspace(char c)
Definition: char.H:53
Input from memory buffer stream.
Definition: IStringStream.H:49
label lineNumber() const
Return current stream line number.
Definition: IOstream.H:438
const dimensionedScalar c
Speed of light in a vacuum.
The geometricSurfacePatch is like patchIdentifier but for surfaces. Holds type, name and index...
messageStream Info
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
static scalar parseNASCoord(const string &s)
Definition: readNAS.C:50
Namespace for OpenFOAM.
ISstream & getLine(string &)
Raw, low-level getline into a string function.
Definition: ISstreamI.H:77