vtkUnstructuredReader.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "vtkUnstructuredReader.H"
30 #include "labelIOField.H"
31 #include "scalarIOField.H"
32 #include "stringIOList.H"
33 #include "cellModel.H"
34 #include "vectorIOField.H"
35 #include "triPointRef.H"
36 
37 /* * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * */
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(vtkUnstructuredReader, 1);
42 }
43 
44 const Foam::Enum
45 <
47 >
49 ({
50  { vtkDataType::VTK_INT, "int" },
51  { vtkDataType::VTK_UINT, "unsigned_int" },
52  { vtkDataType::VTK_LONG, "long" },
53  { vtkDataType::VTK_ULONG, "unsigned_long" },
54  { vtkDataType::VTK_FLOAT, "float" },
55  { vtkDataType::VTK_DOUBLE, "double" },
56  { vtkDataType::VTK_STRING, "string" },
57  { vtkDataType::VTK_ID, "vtkIdType" },
58 });
59 
60 
61 const Foam::Enum
62 <
64 >
66 ({
67  { vtkDataSetType::VTK_FIELD, "FIELD" },
68  { vtkDataSetType::VTK_SCALARS, "SCALARS" },
69  { vtkDataSetType::VTK_VECTORS, "VECTORS" },
70 });
71 
72 
73 const Foam::Enum
74 <
76 >
78 ({
79  { parseMode::NOMODE, "NOMODE" },
80  { parseMode::UNSTRUCTURED_GRID, "UNSTRUCTURED_GRID" },
81  { parseMode::POLYDATA, "POLYDATA" },
82  { parseMode::CELL_DATA, "CELL_DATA" },
83  { parseMode::POINT_DATA, "POINT_DATA" },
84 });
85 
86 
87 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
88 
89 void Foam::vtkUnstructuredReader::warnUnhandledType
90 (
91  Istream& inFile,
92  const label type,
93  labelHashSet& warningGiven
94 ) const
95 {
96  if (warningGiven.insert(type))
97  {
98  IOWarningInFunction(inFile)
99  << "Skipping unknown cell type " << type << nl;
100  }
101 }
102 
103 
104 void Foam::vtkUnstructuredReader::extractCells
105 (
106  Istream& inFile,
107  const labelList& cellTypes,
108  const labelList& cellVertData
109 )
110 {
111  const cellModel& hex = cellModel::ref(cellModel::HEX);
112  const cellModel& prism = cellModel::ref(cellModel::PRISM);
113  const cellModel& pyr = cellModel::ref(cellModel::PYR);
114  const cellModel& tet = cellModel::ref(cellModel::TET);
115 
116  labelList tetPoints(4);
117  labelList pyrPoints(5);
118  labelList prismPoints(6);
119  labelList hexPoints(8);
120 
121  label celli = cells_.size();
122  cells_.setSize(celli+cellTypes.size());
123  cellMap_.setSize(cells_.size(), -1);
124 
125  label facei = faces_.size();
126  faces_.setSize(facei+cellTypes.size());
127  faceMap_.setSize(faces_.size(), -1);
128 
129  label lineI = lines_.size();
130  lines_.setSize(lineI+cellTypes.size());
131  lineMap_.setSize(lines_.size(), -1);
132 
133  label dataIndex = 0;
134 
135 
136  // To mark whether unhandled type has been visited.
137  labelHashSet warningGiven;
138 
139  forAll(cellTypes, i)
140  {
141  switch (cellTypes[i])
142  {
144  {
145  warnUnhandledType(inFile, cellTypes[i], warningGiven);
146  label nRead = cellVertData[dataIndex++];
147  if (nRead != 1)
148  {
149  FatalIOErrorInFunction(inFile)
150  << "Expected size 1 for VTK_VERTEX but found "
151  << nRead << exit(FatalIOError);
152  }
153  dataIndex += nRead;
154  }
155  break;
156 
158  {
159  warnUnhandledType(inFile, cellTypes[i], warningGiven);
160  label nRead = cellVertData[dataIndex++];
161  dataIndex += nRead;
162  }
163  break;
164 
166  {
167  //warnUnhandledType(inFile, cellTypes[i], warningGiven);
168  label nRead = cellVertData[dataIndex++];
169  if (nRead != 2)
170  {
171  FatalIOErrorInFunction(inFile)
172  << "Expected size 2 for VTK_LINE but found "
173  << nRead << exit(FatalIOError);
174  }
175  lineMap_[lineI] = i;
176  labelList& segment = lines_[lineI++];
177  segment.setSize(2);
178  segment[0] = cellVertData[dataIndex++];
179  segment[1] = cellVertData[dataIndex++];
180  }
181  break;
182 
184  {
185  //warnUnhandledType(inFile, cellTypes[i], warningGiven);
186  label nRead = cellVertData[dataIndex++];
187  lineMap_[lineI] = i;
188  labelList& segment = lines_[lineI++];
189  segment.setSize(nRead);
190  for (label& pointi : segment)
191  {
192  pointi = cellVertData[dataIndex++];
193  }
194  }
195  break;
196 
198  {
199  faceMap_[facei] = i;
200  face& f = faces_[facei++];
201  f.setSize(3);
202  label nRead = cellVertData[dataIndex++];
203  if (nRead != 3)
204  {
205  FatalIOErrorInFunction(inFile)
206  << "Expected size 3 for VTK_TRIANGLE but found "
207  << nRead << exit(FatalIOError);
208  }
209  f[0] = cellVertData[dataIndex++];
210  f[1] = cellVertData[dataIndex++];
211  f[2] = cellVertData[dataIndex++];
212  }
213  break;
214 
216  {
217  faceMap_[facei] = i;
218  face& f = faces_[facei++];
219  f.setSize(4);
220  label nRead = cellVertData[dataIndex++];
221  if (nRead != 4)
222  {
223  FatalIOErrorInFunction(inFile)
224  << "Expected size 4 for VTK_QUAD but found "
225  << nRead << exit(FatalIOError);
226  }
227  f[0] = cellVertData[dataIndex++];
228  f[1] = cellVertData[dataIndex++];
229  f[2] = cellVertData[dataIndex++];
230  f[3] = cellVertData[dataIndex++];
231  }
232  break;
233 
235  {
236  faceMap_[facei] = i;
237  face& f = faces_[facei++];
238  label nRead = cellVertData[dataIndex++];
239  f.setSize(nRead);
240  for (label& pointi : f)
241  {
242  pointi = cellVertData[dataIndex++];
243  }
244  }
245  break;
246 
248  {
249  label nRead = cellVertData[dataIndex++];
250  if (nRead != 4)
251  {
252  FatalIOErrorInFunction(inFile)
253  << "Expected size 4 for VTK_TETRA but found "
254  << nRead << exit(FatalIOError);
255  }
256  tetPoints[0] = cellVertData[dataIndex++];
257  tetPoints[1] = cellVertData[dataIndex++];
258  tetPoints[2] = cellVertData[dataIndex++];
259  tetPoints[3] = cellVertData[dataIndex++];
260  cellMap_[celli] = i;
261  cells_[celli++] = cellShape(tet, tetPoints, true);
262  }
263  break;
264 
266  {
267  label nRead = cellVertData[dataIndex++];
268  if (nRead != 5)
269  {
270  FatalIOErrorInFunction(inFile)
271  << "Expected size 5 for VTK_PYRAMID but found "
272  << nRead << exit(FatalIOError);
273  }
274  pyrPoints[0] = cellVertData[dataIndex++];
275  pyrPoints[1] = cellVertData[dataIndex++];
276  pyrPoints[2] = cellVertData[dataIndex++];
277  pyrPoints[3] = cellVertData[dataIndex++];
278  pyrPoints[4] = cellVertData[dataIndex++];
279  cellMap_[celli] = i;
280  cells_[celli++] = cellShape(pyr, pyrPoints, true);
281  }
282  break;
283 
285  {
286  label nRead = cellVertData[dataIndex++];
287  if (nRead != 6)
288  {
289  FatalIOErrorInFunction(inFile)
290  << "Expected size 6 for VTK_WEDGE but found "
291  << nRead << exit(FatalIOError);
292  }
293  //- From mesh description in vtk documentation
294  prismPoints[0] = cellVertData[dataIndex++];
295  prismPoints[2] = cellVertData[dataIndex++];
296  prismPoints[1] = cellVertData[dataIndex++];
297  prismPoints[3] = cellVertData[dataIndex++];
298  prismPoints[5] = cellVertData[dataIndex++];
299  prismPoints[4] = cellVertData[dataIndex++];
300  cellMap_[celli] = i;
301  cells_[celli++] = cellShape(prism, prismPoints, true);
302  }
303  break;
304 
306  {
307  label nRead = cellVertData[dataIndex++];
308  if (nRead != 8)
309  {
310  FatalIOErrorInFunction(inFile)
311  << "Expected size 8 for VTK_HEXAHEDRON but found "
312  << nRead << exit(FatalIOError);
313  }
314  hexPoints[0] = cellVertData[dataIndex++];
315  hexPoints[1] = cellVertData[dataIndex++];
316  hexPoints[2] = cellVertData[dataIndex++];
317  hexPoints[3] = cellVertData[dataIndex++];
318  hexPoints[4] = cellVertData[dataIndex++];
319  hexPoints[5] = cellVertData[dataIndex++];
320  hexPoints[6] = cellVertData[dataIndex++];
321  hexPoints[7] = cellVertData[dataIndex++];
322  cellMap_[celli] = i;
323  cells_[celli++] = cellShape(hex, hexPoints, true);
324  }
325  break;
326 
327  default:
328  warnUnhandledType(inFile, cellTypes[i], warningGiven);
329  label nRead = cellVertData[dataIndex++];
330  dataIndex += nRead;
331  }
332  }
333 
334  DebugInfo
335  << "Read " << celli << " cells;" << facei << " faces." << nl;
336 
337  cells_.setSize(celli);
338  cellMap_.setSize(celli);
339  faces_.setSize(facei);
340  faceMap_.setSize(facei);
341  lines_.setSize(lineI);
342  lineMap_.setSize(lineI);
343 }
344 
345 
346 void Foam::vtkUnstructuredReader::readField
347 (
348  ISstream& inFile,
349  objectRegistry& obj,
350  const word& arrayName,
351  const word& dataType,
352  const label size
353 ) const
354 {
355  if (vtkDataTypeNames.found(dataType))
356  {
357  switch (vtkDataTypeNames[dataType])
358  {
359  case VTK_INT:
360  case VTK_UINT:
361  case VTK_LONG:
362  case VTK_ULONG:
363  case VTK_ID:
364  {
365  auto fieldVals = autoPtr<labelIOField>::New
366  (
367  IOobject(arrayName, "", obj),
368  size
369  );
370  readBlock(inFile, fieldVals().size(), fieldVals());
371  regIOobject::store(fieldVals);
372  }
373  break;
374 
375  case VTK_FLOAT:
376  case VTK_DOUBLE:
377  {
378  auto fieldVals = autoPtr<scalarIOField>::New
379  (
380  IOobject(arrayName, "", obj),
381  size
382  );
383  readBlock(inFile, fieldVals().size(), fieldVals());
384  regIOobject::store(fieldVals);
385  }
386  break;
387 
388  case VTK_STRING:
389  {
390  DebugInfo
391  << "Reading strings:" << size << nl;
392 
393  auto fieldVals = autoPtr<stringIOList>::New
394  (
395  IOobject(arrayName, "", obj),
396  size
397  );
398  // Consume current line.
399  inFile.getLine(fieldVals()[0]);
400 
401  // Read without parsing
402  for (string& s : fieldVals())
403  {
404  inFile.getLine(s);
405  }
406  regIOobject::store(fieldVals);
407  }
408  break;
409 
410  default:
411  {
412  IOWarningInFunction(inFile)
413  << "Unhandled type " << dataType << nl
414  << "Skipping " << size
415  << " words." << nl;
416  scalarField fieldVals;
417  readBlock(inFile, size, fieldVals);
418  }
419  break;
420  }
421  }
422  else
423  {
424  IOWarningInFunction(inFile)
425  << "Unhandled type " << dataType << nl
426  << "Skipping " << size
427  << " words." << nl;
428  scalarField fieldVals;
429  readBlock(inFile, size, fieldVals);
430  }
431 }
432 
433 
434 Foam::wordList Foam::vtkUnstructuredReader::readFieldArray
435 (
436  ISstream& inFile,
437  objectRegistry& obj,
438  const label wantedSize
439 ) const
440 {
441  DynamicList<word> fields;
442 
443  word dataName(inFile);
444 
445  DebugInfo
446  << "dataName:" << dataName << nl;
447 
448  label numArrays(readLabel(inFile));
449  if (debug)
450  {
451  Pout<< "numArrays:" << numArrays << nl;
452  }
453  for (label i = 0; i < numArrays; i++)
454  {
455  word arrayName(inFile);
456  label numComp(readLabel(inFile));
457  label numTuples(readLabel(inFile));
458  word dataType(inFile);
459 
460  DebugInfo
461  << "Reading field " << arrayName
462  << " of " << numTuples << " tuples of rank " << numComp << nl;
463 
464  if (wantedSize != -1 && numTuples != wantedSize)
465  {
466  FatalIOErrorInFunction(inFile)
467  << "Expected " << wantedSize << " tuples but only have "
468  << numTuples << exit(FatalIOError);
469  }
470 
471  readField
472  (
473  inFile,
474  obj,
475  arrayName,
476  dataType,
477  numTuples*numComp
478  );
479  fields.append(arrayName);
480  }
481  return fields.shrink();
482 }
483 
484 
485 Foam::objectRegistry& Foam::vtkUnstructuredReader::selectRegistry
486 (
487  const parseMode readMode
488 )
489 {
490  if (readMode == CELL_DATA)
491  {
492  return cellData_;
493  }
494  else if (readMode == POINT_DATA)
495  {
496  return pointData_;
497  }
498 
499  return otherData_;
500 }
501 
502 
503 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
504 
506 (
507  const objectRegistry& obr,
508  ISstream& inFile
509 )
510 :
511  cellData_(IOobject("cellData", obr)),
512  pointData_(IOobject("pointData", obr)),
513  otherData_(IOobject("otherData", obr))
514 {
515  read(inFile);
516 }
517 
518 
519 void Foam::vtkUnstructuredReader::read(ISstream& inFile)
520 {
521  inFile.getLine(header_);
522  DebugInfo<< "Header : " << header_ << nl;
523 
524  inFile.getLine(title_);
525  DebugInfo<< "Title : " << title_ << nl;
526 
527  inFile.getLine(dataType_);
528  DebugInfo<< "dataType : " << dataType_ << nl;
529 
530  if (dataType_ == "BINARY")
531  {
532  FatalIOErrorInFunction(inFile)
533  << "Binary reading not supported"
534  << exit(FatalIOError);
535  }
536 
537  parseMode readMode = NOMODE;
538  label wantedSize = -1;
539 
540 
541  // Temporary storage for vertices of cells.
542  labelList cellVerts;
543 
544  token tok;
545 
546  while (inFile.read(tok).good() && tok.isWord())
547  {
548  const word tag = tok.wordToken();
549 
550  DebugInfo
551  << "line:" << inFile.lineNumber()
552  << " tag:" << tag << nl;
553 
554  if (tag == "DATASET")
555  {
556  word geomType(inFile);
557  DebugInfo<< "geomType : " << geomType << nl;
558 
559  readMode = parseModeNames[geomType];
560  wantedSize = -1;
561  }
562  else if (tag == "POINTS")
563  {
564  label nPoints(readLabel(inFile));
565  points_.setSize(nPoints);
566 
567  DebugInfo
568  << "Reading " << nPoints << " numbers representing "
569  << points_.size() << " coordinates." << nl;
570 
571  word primitiveTag(inFile);
572  if (primitiveTag != "float" && primitiveTag != "double")
573  {
574  FatalIOErrorInFunction(inFile)
575  << "Expected 'float' entry but found "
576  << primitiveTag
577  << exit(FatalIOError);
578  }
579  for (point& p : points_)
580  {
581  inFile >> p.x() >> p.y() >> p.z();
582  }
583  }
584  else if (tag == "CELLS")
585  {
586  label nCells(readLabel(inFile));
587  label nNumbers(readLabel(inFile));
588  DebugInfo
589  << "Reading " << nCells << " cells or faces." << nl;
590 
591  readBlock(inFile, nNumbers, cellVerts);
592  }
593  else if (tag == "CELL_TYPES")
594  {
595  label nCellTypes(readLabel(inFile));
596 
598  readBlock(inFile, nCellTypes, cellTypes);
599 
600  if (cellTypes.size() > 0 && cellVerts.size() == 0)
601  {
602  FatalIOErrorInFunction(inFile)
603  << "Found " << cellTypes.size()
604  << " cellTypes but no cells."
605  << exit(FatalIOError);
606  }
607 
608  extractCells(inFile, cellTypes, cellVerts);
609  cellVerts.clear();
610  }
611  else if (tag == "LINES")
612  {
613  label nLines(readLabel(inFile));
614  label nNumbers(readLabel(inFile));
615  DebugInfo
616  << "Reading " << nLines << " lines." << nl;
617 
618  labelList lineVerts;
619  readBlock(inFile, nNumbers, lineVerts);
620 
621  label lineI = lines_.size();
622  lines_.setSize(lineI+nLines);
623  lineMap_.setSize(lines_.size());
624 
625  label elemI = 0;
626  for (label i = 0; i < nLines; i++)
627  {
628  lineMap_[lineI] = lineI;
629  labelList& f = lines_[lineI];
630  f.setSize(lineVerts[elemI++]);
631  for (label& pointi : f)
632  {
633  pointi = lineVerts[elemI++];
634  }
635  lineI++;
636  }
637  }
638  else if (tag == "POLYGONS")
639  {
640  // If in polydata mode
641 
642  label nFaces(readLabel(inFile));
643  label nNumbers(readLabel(inFile));
644  DebugInfo
645  << "Reading " << nFaces << " faces." << nl;
646 
647  labelList faceVerts;
648  readBlock(inFile, nNumbers, faceVerts);
649 
650  label facei = faces_.size();
651  faces_.setSize(facei+nFaces);
652  faceMap_.setSize(faces_.size());
653 
654  label elemI = 0;
655  for (label i = 0; i < nFaces; i++)
656  {
657  faceMap_[facei] = facei;
658  face& f = faces_[facei];
659  f.setSize(faceVerts[elemI++]);
660  for (label& pointi : f)
661  {
662  pointi = faceVerts[elemI++];
663  }
664  facei++;
665  }
666  }
667  else if (tag == "POINT_DATA")
668  {
669  // 'POINT_DATA 24'
670  readMode = POINT_DATA;
671  wantedSize = points_.size();
672 
673  label nPoints(readLabel(inFile));
674  if (nPoints != wantedSize)
675  {
676  FatalIOErrorInFunction(inFile)
677  << "Reading POINT_DATA : expected " << wantedSize
678  << " but read " << nPoints << exit(FatalIOError);
679  }
680  }
681  else if (tag == "CELL_DATA")
682  {
683  readMode = CELL_DATA;
684  wantedSize = cells_.size()+faces_.size()+lines_.size();
685 
686  label nCells(readLabel(inFile));
687  if (nCells != wantedSize)
688  {
689  FatalIOErrorInFunction(inFile)
690  << "Reading CELL_DATA : expected "
691  << wantedSize
692  << " but read " << nCells << exit(FatalIOError);
693  }
694  }
695  else if (tag == "FIELD")
696  {
697  // wantedSize already set according to type we expected to read.
698  readFieldArray(inFile, selectRegistry(readMode), wantedSize);
699  }
700  else if (tag == "SCALARS")
701  {
702  string line;
703  inFile.getLine(line);
704  IStringStream is(line);
705  word dataName(is);
706  word dataType(is);
707  //label numComp(readLabel(inFile));
708 
709  DebugInfo
710  << "Reading scalar " << dataName
711  << " of type " << dataType
712  << " from lookup table" << nl;
713 
714  word lookupTableTag(inFile);
715  if (lookupTableTag != "LOOKUP_TABLE")
716  {
717  FatalIOErrorInFunction(inFile)
718  << "Expected tag LOOKUP_TABLE but read "
719  << lookupTableTag
720  << exit(FatalIOError);
721  }
722 
723  word lookupTableName(inFile);
724 
725  readField
726  (
727  inFile,
728  selectRegistry(readMode),
729  dataName,
730  dataType,
731  wantedSize//*numComp
732  );
733  }
734  else if (tag == "VECTORS" || tag == "NORMALS")
735  {
736  // 'NORMALS Normals float'
737  string line;
738  inFile.getLine(line);
739  IStringStream is(line);
740  word dataName(is);
741  word dataType(is);
742  DebugInfo
743  << "Reading vector " << dataName
744  << " of type " << dataType << nl;
745 
746  objectRegistry& reg = selectRegistry(readMode);
747 
748  readField
749  (
750  inFile,
751  reg,
752  dataName,
753  dataType,
754  3*wantedSize
755  );
756 
757  if
758  (
759  vtkDataTypeNames[dataType] == VTK_FLOAT
760  || vtkDataTypeNames[dataType] == VTK_DOUBLE
761  )
762  {
763  objectRegistry::iterator iter = reg.find(dataName);
764  scalarField s(*dynamic_cast<const scalarField*>(iter()));
765  reg.erase(iter);
766  auto fieldVals = autoPtr<vectorIOField>::New
767  (
768  IOobject(dataName, "", reg),
769  s.size()/3
770  );
771 
772  label elemI = 0;
773  for (vector& val : fieldVals())
774  {
775  val.x() = s[elemI++];
776  val.y() = s[elemI++];
777  val.z() = s[elemI++];
778  }
779  regIOobject::store(fieldVals);
780  }
781  }
782  else if (tag == "TEXTURE_COORDINATES")
783  {
784  // 'TEXTURE_COORDINATES TCoords 2 float'
785  string line;
786  inFile.getLine(line);
787  IStringStream is(line);
788  word dataName(is); //"Tcoords"
789  label dim(readLabel(is));
790  word dataType(is);
791 
792  DebugInfo
793  << "Reading texture coords " << dataName
794  << " dimension " << dim
795  << " of type " << dataType << nl;
796 
797  scalarField coords(dim*points_.size());
798  readBlock(inFile, coords.size(), coords);
799  }
800  else if (tag == "TRIANGLE_STRIPS")
801  {
802  label nStrips(readLabel(inFile));
803  label nNumbers(readLabel(inFile));
804  DebugInfo
805  << "Reading " << nStrips << " triangle strips." << nl;
806 
807  labelList faceVerts;
808  readBlock(inFile, nNumbers, faceVerts);
809 
810  // Count number of triangles
811  label elemI = 0;
812  label nTris = 0;
813  for (label i = 0; i < nStrips; i++)
814  {
815  label nVerts = faceVerts[elemI++];
816  nTris += nVerts-2;
817  elemI += nVerts;
818  }
819 
820 
821  // Store
822  label facei = faces_.size();
823  faces_.setSize(facei+nTris);
824  faceMap_.setSize(faces_.size());
825  elemI = 0;
826  for (label i = 0; i < nStrips; i++)
827  {
828  label nVerts = faceVerts[elemI++];
829  label nTris = nVerts-2;
830 
831  // Read first triangle
832  faceMap_[facei] = facei;
833  face& f = faces_[facei++];
834  f.setSize(3);
835  f[0] = faceVerts[elemI++];
836  f[1] = faceVerts[elemI++];
837  f[2] = faceVerts[elemI++];
838  for (label triI = 1; triI < nTris; triI++)
839  {
840  faceMap_[facei] = facei;
841  face& f = faces_[facei++];
842  f.setSize(3);
843  f[0] = faceVerts[elemI-1];
844  f[1] = faceVerts[elemI-2];
845  f[2] = faceVerts[elemI++];
846  }
847  }
848  }
849  else if (tag == "METADATA")
850  {
851  word infoTag(inFile);
852  if (infoTag != "INFORMATION")
853  {
854  FatalIOErrorInFunction(inFile)
855  << "Unsupported tag "
856  << infoTag << exit(FatalIOError);
857  }
858  label nInfo(readLabel(inFile));
859  DebugInfo
860  << "Consuming " << nInfo << " metadata information." << nl;
861 
862  string line;
863  // Consume rest of line
864  inFile.getLine(line);
865  for (label i = 0; i < 2*nInfo; i++)
866  {
867  inFile.getLine(line);
868  }
869  }
870  else
871  {
872  FatalIOErrorInFunction(inFile)
873  << "Unsupported tag "
874  << tag << exit(FatalIOError);
875  }
876  }
877 
878 
879  // There is some problem with orientation of prisms - the point
880  // ordering seems to be different for some exports (e.g. of cgns)
881  {
882  const cellModel& prism = cellModel::ref(cellModel::PRISM);
883 
884  label nSwapped = 0;
885 
886  for (cellShape& shape : cells_)
887  {
888  if (shape.model() == prism)
889  {
890  const triPointRef bottom
891  (
892  points_[shape[0]],
893  points_[shape[1]],
894  points_[shape[2]]
895  );
896  const triPointRef top
897  (
898  points_[shape[3]],
899  points_[shape[4]],
900  points_[shape[5]]
901  );
902 
903  const point bottomCc(bottom.centre());
904  const vector bottomNormal(bottom.areaNormal());
905  const point topCc(top.centre());
906 
907  if (((topCc - bottomCc) & bottomNormal) < 0)
908  {
909  // Flip top and bottom
910  Swap(shape[0], shape[3]);
911  Swap(shape[1], shape[4]);
912  Swap(shape[2], shape[5]);
913  nSwapped++;
914  }
915  }
916  }
917  if (nSwapped > 0)
918  {
919  WarningInFunction << "Swapped " << nSwapped << " prismatic cells"
920  << nl;
921  }
922  }
923 
924  if (debug)
925  {
926  Info<< "Read points:" << points_.size()
927  << " cellShapes:" << cells_.size()
928  << " faces:" << faces_.size()
929  << " lines:" << lines_.size()
930  << nl << nl;
931 
932  Info<< "Cell fields:" << nl;
933  printFieldStats<vectorIOField>(cellData_);
934  printFieldStats<scalarIOField>(cellData_);
935  printFieldStats<labelIOField>(cellData_);
936  printFieldStats<stringIOList>(cellData_);
937  Info<< nl << nl;
938 
939  Info<< "Point fields:" << nl;
940  printFieldStats<vectorIOField>(pointData_);
941  printFieldStats<scalarIOField>(pointData_);
942  printFieldStats<labelIOField>(pointData_);
943  printFieldStats<stringIOList>(pointData_);
944  Info<< nl << nl;
945 
946  Info<< "Other fields:" << nl;
947  printFieldStats<vectorIOField>(otherData_);
948  printFieldStats<scalarIOField>(otherData_);
949  printFieldStats<labelIOField>(otherData_);
950  printFieldStats<stringIOList>(otherData_);
951  }
952 }
953 
954 
955 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::vtkUnstructuredReader::vtkDataSetTypeNames
static const Enum< vtkDataSetType > vtkDataSetTypeNames
Definition: vtkUnstructuredReader.H:96
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::val
label ListType::const_reference val
Definition: ListOps.H:407
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::vtkUnstructuredReader::vtkDataSetType
vtkDataSetType
Enumeration defining the vtk dataset types.
Definition: vtkUnstructuredReader.H:89
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
Foam::ISstream::getLine
ISstream & getLine(std::string &str, char delim='\n')
Raw, low-level getline into a string function.
Definition: ISstreamI.H:78
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::vtkUnstructuredReader::parseModeNames
static const Enum< parseMode > parseModeNames
Definition: vtkUnstructuredReader.H:109
vtkUnstructuredReader.H
Foam::expressions::patchExpr::POINT_DATA
Point data.
Definition: patchExprFwd.H:60
Foam::cellModel::HEX
hex
Definition: cellModel.H:81
Foam::vtk::VTK_HEXAHEDRON
Definition: foamVtkCore.H:103
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::ISstream
Generic input stream using standard (STL) streams.
Definition: ISstream.H:54
Foam::vtk::VTK_PYRAMID
Definition: foamVtkCore.H:105
Foam::FatalIOError
IOerror FatalIOError
scalarIOField.H
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
triPointRef.H
Foam::Swap
void Swap(DynamicList< T, SizeMin1 > &a, DynamicList< T, SizeMin2 > &b)
Definition: DynamicListI.H:909
Foam::cellModel::TET
tet
Definition: cellModel.H:85
Foam::vtk::VTK_LINE
Definition: foamVtkCore.H:94
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::vtkUnstructuredReader::vtkUnstructuredReader
vtkUnstructuredReader(const objectRegistry &obr, ISstream &)
Construct from Istream, read all.
Definition: vtkUnstructuredReader.C:506
Foam::vtkUnstructuredReader::vtkDataType
vtkDataType
Enumeration defining the vtk data types.
Definition: vtkUnstructuredReader.H:73
Foam::label
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:62
cellModel.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::vtk::VTK_POLY_VERTEX
Definition: foamVtkCore.H:93
Foam::vtk::VTK_WEDGE
Definition: foamVtkCore.H:104
Foam::cellModel::PRISM
prism
Definition: cellModel.H:83
Foam::cellModel::ref
static const cellModel & ref(const modelType model)
Look up reference to cellModel by enumeration. Fatal on failure.
Definition: cellModels.C:157
Foam::vtk::VTK_TRIANGLE
Definition: foamVtkCore.H:96
Foam::segment
Pair< point > segment
Definition: distributedTriSurfaceMesh.H:75
Foam::vtk::VTK_POLY_LINE
Definition: foamVtkCore.H:95
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vtk::VTK_VERTEX
Definition: foamVtkCore.H:92
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::hex
IOstream & hex(IOstream &io)
Definition: IOstream.H:429
vectorIOField.H
Foam::regIOobject::store
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:37
Foam::cellModel::PYR
pyr
Definition: cellModel.H:84
cellTypes
const labelList & cellTypes
Definition: setCellMask.H:33
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::nl
constexpr char nl
Definition: Ostream.H:372
f
labelList f(nPoints)
Foam::readLabel
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition: label.H:70
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:102
Foam::vtk::VTK_POLYGON
Definition: foamVtkCore.H:98
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::vtk::VTK_TETRA
Definition: foamVtkCore.H:101
Foam::ISstream::read
virtual Istream & read(token &t)
Return next token from stream.
Definition: ISstream.C:158
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::IOstream::lineNumber
label lineNumber() const
Const access to the current stream line number.
Definition: IOstream.H:301
Foam::vtk::fileTag::CELL_DATA
"CellData"
Foam::vtk::VTK_QUAD
Definition: foamVtkCore.H:100
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
IOWarningInFunction
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Definition: messageStream.H:306
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:216
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:415
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::vtkUnstructuredReader::vtkDataTypeNames
static const Enum< vtkDataType > vtkDataTypeNames
Definition: vtkUnstructuredReader.H:85
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::vtkUnstructuredReader::parseMode
parseMode
Enumeration defining the parse mode - type of data being read.
Definition: vtkUnstructuredReader.H:100
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:78
labelIOField.H
stringIOList.H