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-2024 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 "triangle.H"
36 #include "stringOps.H"
37 
38 /* * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * */
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(vtkUnstructuredReader, 1);
43 }
44 
45 const Foam::Enum
46 <
48 >
50 ({
51  { vtkDataType::VTK_INT, "int" },
52  // Not yet required: { vtkDataType::VTK_INT64, "vtktypeint64" },
53  { vtkDataType::VTK_UINT, "unsigned_int" },
54  { vtkDataType::VTK_LONG, "long" },
55  { vtkDataType::VTK_ULONG, "unsigned_long" },
56  { vtkDataType::VTK_FLOAT, "float" },
57  { vtkDataType::VTK_DOUBLE, "double" },
58  { vtkDataType::VTK_STRING, "string" },
59  { vtkDataType::VTK_ID, "vtkIdType" },
60 });
61 
62 
63 const Foam::Enum
64 <
66 >
68 ({
69  { vtkDataSetType::VTK_FIELD, "FIELD" },
70  { vtkDataSetType::VTK_SCALARS, "SCALARS" },
71  { vtkDataSetType::VTK_VECTORS, "VECTORS" },
72 });
73 
74 
75 const Foam::Enum
76 <
78 >
80 ({
81  { parseMode::NOMODE, "NOMODE" },
82  { parseMode::UNSTRUCTURED_GRID, "UNSTRUCTURED_GRID" },
83  { parseMode::POLYDATA, "POLYDATA" },
84  { parseMode::CELL_DATA, "CELL_DATA" },
85  { parseMode::POINT_DATA, "POINT_DATA" },
86 });
87 
88 
89 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
90 
91 namespace Foam
92 {
93 
94 // Read N elements into List
95 template<class T>
96 static inline void readBlock(Istream& is, const label n, List<T>& list)
97 {
98  list.resize(n);
99  for (T& val : list)
100  {
101  is >> val;
102  }
103 }
104 
105 } // End namespace Foam
106 
107 
108 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
109 
110 void Foam::vtkUnstructuredReader::warnUnhandledType
111 (
112  const Istream& is,
113  const label type,
114  labelHashSet& warningGiven
115 )
116 {
117  if (warningGiven.insert(type))
118  {
120  << "Skipping unknown cell type " << type << nl;
121  }
122 }
123 
124 
125 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
126 
127 void Foam::vtkUnstructuredReader::readOffsetsConnectivity
128 (
129  ISstream& is,
130  const char* entryName,
131  const label nOffsets,
132  labelList& offsets,
133  const label nConnectivity,
134  labelList& connectivity
135 )
136 {
137  token tok;
138 
139  is.read(tok);
140  if (!tok.isWord("OFFSETS"))
141  {
143  << "Expected OFFSETS for " << entryName
144  << ", found "
145  << tok.info() << nl
146  << exit(FatalIOError);
147  }
148  is.getLine(nullptr); // Consume rest of line
149  readBlock(is, nOffsets, offsets);
150 
151  is.read(tok);
152  if (!tok.isWord("CONNECTIVITY"))
153  {
155  << "Expected CONNECTIVITY for " << entryName
156  << ", found " << tok.info() << nl
157  << exit(FatalIOError);
158  }
159  is.getLine(nullptr); // Consume rest of line
160  readBlock(is, nConnectivity, connectivity);
161 }
162 
163 
164 void Foam::vtkUnstructuredReader::extractCells
165 (
166  const Istream& is,
167  const labelUList& cellTypes,
168  const labelUList& elemOffsets,
169  const labelUList& elemVerts
170 )
171 {
172  const cellModel& hex = cellModel::ref(cellModel::HEX);
173  const cellModel& prism = cellModel::ref(cellModel::PRISM);
174  const cellModel& pyr = cellModel::ref(cellModel::PYR);
175  const cellModel& tet = cellModel::ref(cellModel::TET);
176 
177  // Pass 0: sizing
178  label nCells = 0, nFaces = 0, nLines = 0;
179  forAll(cellTypes, elemi)
180  {
181  switch (cellTypes[elemi])
182  {
185  {
186  ++nLines;
187  }
188  break;
189 
193  {
194  ++nFaces;
195  }
196  break;
197 
202  {
203  ++nCells;
204  }
205  break;
206 
207  default:
208  break;
209  }
210  }
211 
212  label celli = cells_.size();
213  cells_.resize(celli + nCells);
214  cellMap_.resize(cells_.size(), -1);
215 
216  label facei = faces_.size();
217  faces_.resize(facei + nFaces);
218  faceMap_.resize(faces_.size(), -1);
219 
220  label linei = lines_.size();
221  lines_.resize(linei + nLines);
222  lineMap_.resize(lines_.size(), -1);
223 
224  // General scratch space for cell vertices
225  labelList cellPoints(16);
226 
227  label dataIndex = 0; // Addressing into vertices stream
228 
229  // To mark whether unhandled type has been visited.
230  labelHashSet warningGiven;
231 
232 
233  forAll(cellTypes, elemi)
234  {
235  // Vertices per element - from offsets or embedded size
236  const label nVerts =
237  (
238  elemOffsets.empty()
239  ? elemVerts[dataIndex++]
240  : (elemOffsets[elemi+1] - elemOffsets[elemi])
241  );
242 
243  // The addresseed vertices
244  const labelSubList verts(elemVerts, nVerts, dataIndex);
245  dataIndex += nVerts;
246 
247  switch (cellTypes[elemi])
248  {
250  {
251  warnUnhandledType(is, cellTypes[elemi], warningGiven);
252  if (nVerts != 1)
253  {
255  << "Expected size 1 for VTK_VERTEX, found "
256  << nVerts << nl
257  << exit(FatalIOError);
258  }
259  }
260  break;
261 
263  {
264  warnUnhandledType(is, cellTypes[elemi], warningGiven);
265  }
266  break;
267 
269  {
270  if (nVerts != 2)
271  {
273  << "Expected size 2 for VTK_LINE, found "
274  << nVerts << nl
275  << exit(FatalIOError);
276  }
277  lineMap_[linei] = elemi;
278 
279  // Same vertex ordering
280  lines_[linei++] = verts;
281  }
282  break;
283 
285  {
286  lineMap_[linei] = elemi;
287 
288  // Same vertex ordering
289  lines_[linei++] = verts;
290  }
291  break;
292 
294  {
295  if (nVerts != 3)
296  {
298  << "Expected size 3 for VTK_TRIANGLE, found "
299  << nVerts << nl
300  << exit(FatalIOError);
301  }
302  faceMap_[facei] = elemi;
303 
304  // Same vertex ordering
305  static_cast<labelList&>(faces_[facei++]) = verts;
306  }
307  break;
308 
310  {
311  if (nVerts != 4)
312  {
314  << "Expected size 4 for VTK_QUAD, found "
315  << nVerts << nl
316  << exit(FatalIOError);
317  }
318  faceMap_[facei] = elemi;
319 
320  // Same vertex ordering
321  static_cast<labelList&>(faces_[facei++]) = verts;
322  }
323  break;
324 
326  {
327  faceMap_[facei] = elemi;
328 
329  // Same vertex ordering
330  static_cast<labelList&>(faces_[facei++]) = verts;
331  }
332  break;
333 
335  {
336  if (nVerts != 4)
337  {
339  << "Expected size 4 for VTK_TETRA, found "
340  << nVerts << nl
341  << exit(FatalIOError);
342  }
343  cellMap_[celli] = elemi;
344 
345  // Same vertex ordering
346  cells_[celli++].reset(tet, verts, true);
347  }
348  break;
349 
351  {
352  if (nVerts != 5)
353  {
355  << "Expected size 5 for VTK_PYRAMID, found "
356  << nVerts << nl
357  << exit(FatalIOError);
358  }
359  cellMap_[celli] = elemi;
360 
361  // Same vertex ordering
362  cells_[celli++].reset(pyr, verts, true);
363  }
364  break;
365 
367  {
368  if (nVerts != 6)
369  {
371  << "Expected size 6 for VTK_WEDGE, found "
372  << nVerts << nl
373  << exit(FatalIOError);
374  }
375  cellMap_[celli] = elemi;
376 
377  // VTK_WEDGE triangles point outwards (swap 1<->2, 4<->5)
378  labelSubList shape(cellPoints, nVerts);
379  shape[0] = verts[0];
380  shape[2] = verts[1];
381  shape[1] = verts[2];
382  shape[3] = verts[3];
383  shape[5] = verts[4];
384  shape[4] = verts[5];
385 
386  cells_[celli++].reset(prism, shape, true);
387  }
388  break;
389 
391  {
392  if (nVerts != 8)
393  {
395  << "Expected size 8 for VTK_HEXAHEDRON, found "
396  << nVerts << nl
397  << exit(FatalIOError);
398  }
399  cellMap_[celli] = elemi;
400 
401  // Same vertex ordering
402  cells_[celli++].reset(hex, verts, true);
403  }
404  break;
405 
406  default:
407  {
408  warnUnhandledType(is, cellTypes[elemi], warningGiven);
409  }
410  break;
411  }
412  }
413 
414  DebugInfo
415  << "Read"
416  << " cells:" << celli
417  << " faces:" << facei
418  << " lines:" << linei
419  << nl;
420 
421  cells_.resize(celli);
422  cellMap_.resize(celli);
423 
424  faces_.resize(facei);
425  faceMap_.resize(facei);
426 
427  lines_.resize(linei);
428  lineMap_.resize(linei);
429 }
430 
431 
432 void Foam::vtkUnstructuredReader::readField
433 (
434  ISstream& inFile,
435  objectRegistry& obj,
436  const word& arrayName,
437  const word& dataType,
438  const label size
439 ) const
440 {
441  if (vtkDataTypeNames.found(dataType))
442  {
443  switch (vtkDataTypeNames[dataType])
444  {
445  case VTK_INT:
446  case VTK_INT64:
447  case VTK_UINT:
448  case VTK_LONG:
449  case VTK_ULONG:
450  case VTK_ID:
451  {
452  auto fieldVals = autoPtr<labelIOField>::New
453  (
454  IOobject(arrayName, "", obj),
455  size
456  );
457  readBlock(inFile, fieldVals().size(), fieldVals());
458  regIOobject::store(fieldVals);
459  }
460  break;
461 
462  case VTK_FLOAT:
463  case VTK_DOUBLE:
464  {
465  auto fieldVals = autoPtr<scalarIOField>::New
466  (
467  IOobject(arrayName, "", obj),
468  size
469  );
470  readBlock(inFile, fieldVals().size(), fieldVals());
471  regIOobject::store(fieldVals);
472  }
473  break;
474 
475  case VTK_STRING:
476  {
477  DebugInfo
478  << "Reading strings:" << size << nl;
479 
480  auto fieldVals = autoPtr<stringIOList>::New
481  (
482  IOobject(arrayName, "", obj),
483  size
484  );
485  // Consume current line.
486  inFile.getLine(fieldVals()[0]);
487 
488  // Read without parsing
489  for (string& s : fieldVals())
490  {
491  inFile.getLine(s);
492  }
493  regIOobject::store(fieldVals);
494  }
495  break;
496 
497  default:
498  {
499  IOWarningInFunction(inFile)
500  << "Unhandled type " << dataType << nl
501  << "Skipping " << size
502  << " words." << nl;
503  scalarField fieldVals;
504  readBlock(inFile, size, fieldVals);
505  }
506  break;
507  }
508  }
509  else
510  {
511  IOWarningInFunction(inFile)
512  << "Unhandled type " << dataType << nl
513  << "Skipping " << size
514  << " words." << nl;
515  scalarField fieldVals;
516  readBlock(inFile, size, fieldVals);
517  }
518 }
519 
520 
521 Foam::wordList Foam::vtkUnstructuredReader::readFieldArray
522 (
523  ISstream& inFile,
524  objectRegistry& obj,
525  const label wantedSize
526 ) const
527 {
528  DynamicList<word> fields;
529 
530  word dataName(inFile);
531 
532  DebugInfo
533  << "dataName:" << dataName << nl;
534 
535  label numArrays(readLabel(inFile));
536  if (debug)
537  {
538  Pout<< "numArrays:" << numArrays << nl;
539  }
540  for (label i = 0; i < numArrays; i++)
541  {
542  word arrayName(inFile);
543  label numComp(readLabel(inFile));
544  label numTuples(readLabel(inFile));
545  word dataType(inFile);
546 
547  DebugInfo
548  << "Reading field " << arrayName
549  << " of " << numTuples << " tuples of rank " << numComp << nl;
550 
551  if (wantedSize != -1 && numTuples != wantedSize)
552  {
553  FatalIOErrorInFunction(inFile)
554  << "Expected " << wantedSize << " tuples but only have "
555  << numTuples << nl
556  << exit(FatalIOError);
557  }
558 
559  readField
560  (
561  inFile,
562  obj,
563  arrayName,
564  dataType,
565  numTuples*numComp
566  );
567  fields.append(arrayName);
568  }
569 
570  return wordList(std::move(fields));
571 }
572 
573 
574 Foam::objectRegistry& Foam::vtkUnstructuredReader::selectRegistry
575 (
576  const parseMode readMode
577 )
578 {
579  if (readMode == CELL_DATA)
580  {
581  return cellData_;
582  }
583  else if (readMode == POINT_DATA)
584  {
585  return pointData_;
586  }
587 
588  return otherData_;
589 }
590 
591 
592 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
593 
595 (
596  const objectRegistry& obr,
597  ISstream& is
598 )
599 :
600  version_(2.0),
601  cellData_(IOobject("cellData", obr)),
602  pointData_(IOobject("pointData", obr)),
603  otherData_(IOobject("otherData", obr))
604 {
605  read(is);
606 }
607 
608 
609 void Foam::vtkUnstructuredReader::read(ISstream& inFile)
610 {
611  inFile.getLine(header_);
612  DebugInfo<< "Header : " << header_ << nl;
613 
614  // Extract version number from "# vtk DataFile Version 5.1"
615  const auto split = stringOps::splitSpace(header_);
616  if (split.size() >= 4 && split[split.size() - 2] == "Version")
617  {
618  float ver(0);
619  if (readFloat(split[split.size() - 1], ver))
620  {
621  version_ = ver;
622  DebugInfo<< "Version : " << version_ << nl;
623  }
624  }
625 
626  inFile.getLine(title_);
627  DebugInfo<< "Title : " << title_ << nl;
628 
629  inFile.getLine(dataType_);
630  DebugInfo<< "dataType : " << dataType_ << nl;
631 
632  if (dataType_ == "BINARY")
633  {
634  FatalIOErrorInFunction(inFile)
635  << "Binary reading not supported" << nl
636  << exit(FatalIOError);
637  }
638 
639  parseMode readMode = NOMODE;
640  label wantedSize = -1;
641 
642 
643  // Intermediate storage for vertices of cells.
644  labelList cellOffsets, cellVerts;
645 
646  token tok;
647 
648  while (inFile.read(tok).good() && tok.isWord())
649  {
650  const word tag = tok.wordToken();
651 
652  DebugInfo
653  << "line:" << inFile.lineNumber()
654  << " tag:" << tag << nl;
655 
656  if (tag == "DATASET")
657  {
658  word geomType(inFile);
659  DebugInfo<< "geomType : " << geomType << nl;
660 
661  readMode = parseModeNames[geomType];
662  wantedSize = -1;
663  }
664  else if (tag == "POINTS")
665  {
666  const label nPoints(readLabel(inFile));
667  points_.resize(nPoints);
668 
669  DebugInfo
670  << "Reading " << nPoints << " coordinates" << nl;
671 
672  word primitiveTag(inFile);
673  if (primitiveTag != "float" && primitiveTag != "double")
674  {
675  FatalIOErrorInFunction(inFile)
676  << "Expected 'float' entry, found "
677  << primitiveTag << nl
678  << exit(FatalIOError);
679  }
680  for (point& p : points_)
681  {
682  inFile >> p.x() >> p.y() >> p.z();
683  }
684  }
685  else if (tag == "CELLS")
686  {
687  label nCells(readLabel(inFile));
688  const label nNumbers(readLabel(inFile));
689 
690  if (version_ < 5.0)
691  {
692  // VTK 4.2 and earlier (single-block)
693  DebugInfo
694  << "Reading " << nCells
695  << " cells/faces (single block)" << nl;
696 
697  cellOffsets.clear();
698  readBlock(inFile, nNumbers, cellVerts);
699  }
700  else
701  {
702  // VTK 5.0 and later (OFFSETS and CONNECTIVITY)
703 
704  const label nOffsets(nCells);
705  --nCells;
706 
707  DebugInfo
708  << "Reading offsets/connectivity for "
709  << nCells << " cells/faces" << nl;
710 
711  readOffsetsConnectivity
712  (
713  inFile,
714  "CELLS",
715  nOffsets, cellOffsets,
716  nNumbers, cellVerts
717  );
718  }
719  }
720  else if (tag == "CELL_TYPES")
721  {
722  const label nCellTypes(readLabel(inFile));
723 
725  readBlock(inFile, nCellTypes, cellTypes);
726 
727  if (!cellTypes.empty() && cellVerts.empty())
728  {
729  FatalIOErrorInFunction(inFile)
730  << "Found " << cellTypes.size()
731  << " cellTypes but no cells." << nl
732  << exit(FatalIOError);
733  }
734 
735  extractCells(inFile, cellTypes, cellOffsets, cellVerts);
736  cellOffsets.clear();
737  cellVerts.clear();
738  }
739  else if (tag == "LINES")
740  {
741  label nLines(readLabel(inFile));
742  const label nNumbers(readLabel(inFile));
743 
744  labelList elemOffsets, elemVerts;
745 
746  if (version_ < 5.0)
747  {
748  // VTK 4.2 and earlier (single-block)
749  DebugInfo
750  << "Reading " << nLines
751  << " lines (single block)" << nl;
752 
753  readBlock(inFile, nNumbers, elemVerts);
754  }
755  else
756  {
757  // VTK 5.0 and later (OFFSETS and CONNECTIVITY)
758 
759  const label nOffsets(nLines);
760  --nLines;
761 
762  DebugInfo
763  << "Reading offsets/connectivity for "
764  << nLines << " lines" << nl;
765 
766  readOffsetsConnectivity
767  (
768  inFile,
769  "LINES",
770  nOffsets, elemOffsets,
771  nNumbers, elemVerts
772  );
773  }
774 
775  // Append into lines
776  label linei = lines_.size();
777  lines_.resize(linei+nLines);
778  lineMap_.resize(lines_.size());
779 
780  label dataIndex = 0;
781  for (label i = 0; i < nLines; i++)
782  {
783  const label nVerts =
784  (
785  elemOffsets.empty()
786  ? elemVerts[dataIndex++]
787  : (elemOffsets[i+1] - elemOffsets[i])
788  );
789  const labelSubList verts(elemVerts, nVerts, dataIndex);
790  dataIndex += nVerts;
791 
792  lineMap_[linei] = linei;
793  lines_[linei++] = verts;
794  }
795  }
796  else if (tag == "POLYGONS")
797  {
798  label nFaces(readLabel(inFile));
799  const label nNumbers(readLabel(inFile));
800 
801  labelList elemOffsets, elemVerts;
802 
803  if (version_ < 5.0)
804  {
805  // VTK 4.2 and earlier (single-block)
806  DebugInfo
807  << "Reading " << nFaces
808  << " faces (single block)" << nl;
809 
810  readBlock(inFile, nNumbers, elemVerts);
811  }
812  else
813  {
814  // VTK 5.0 and later (OFFSETS and CONNECTIVITY)
815 
816  const label nOffsets(nFaces);
817  --nFaces;
818 
819  DebugInfo
820  << "Reading offsets/connectivity for "
821  << nFaces << " faces" << nl;
822 
823  readOffsetsConnectivity
824  (
825  inFile,
826  "POLYGONS",
827  nOffsets, elemOffsets,
828  nNumbers, elemVerts
829  );
830  }
831 
832  // Append into faces
833  label facei = faces_.size();
834  faces_.resize(facei+nFaces);
835  faceMap_.resize(faces_.size());
836 
837  label dataIndex = 0;
838  for (label i = 0; i < nFaces; ++i)
839  {
840  const label nVerts =
841  (
842  elemOffsets.empty()
843  ? elemVerts[dataIndex++]
844  : (elemOffsets[i+1] - elemOffsets[i])
845  );
846  const labelSubList verts(elemVerts, nVerts, dataIndex);
847  dataIndex += nVerts;
848 
849  faceMap_[facei] = facei;
850  static_cast<labelList&>(faces_[facei++]) = verts;
851  }
852  }
853  else if (tag == "POINT_DATA")
854  {
855  // 'POINT_DATA 24'
856  readMode = POINT_DATA;
857  wantedSize = points_.size();
858 
859  label nPoints(readLabel(inFile));
860  if (nPoints != wantedSize)
861  {
862  FatalIOErrorInFunction(inFile)
863  << "Reading POINT_DATA : expected " << wantedSize
864  << " but read " << nPoints << nl
865  << exit(FatalIOError);
866  }
867  }
868  else if (tag == "CELL_DATA")
869  {
870  readMode = CELL_DATA;
871  wantedSize = cells_.size()+faces_.size()+lines_.size();
872 
873  label nCells(readLabel(inFile));
874  if (nCells != wantedSize)
875  {
876  FatalIOErrorInFunction(inFile)
877  << "Reading CELL_DATA : expected "
878  << wantedSize
879  << " but read " << nCells << nl
880  << exit(FatalIOError);
881  }
882  }
883  else if (tag == "FIELD")
884  {
885  // wantedSize already set according to type we expected to read.
886  readFieldArray(inFile, selectRegistry(readMode), wantedSize);
887  }
888  else if (tag == "SCALARS")
889  {
890  string line;
891  inFile.getLine(line);
892  ISpanStream is(line);
893 
894  word dataName(is);
895  word dataType(is);
896  //label numComp(readLabel(inFile));
897 
898  DebugInfo
899  << "Reading scalar " << dataName
900  << " of type " << dataType
901  << " from lookup table" << nl;
902 
903  word lookupTableTag(inFile);
904  if (lookupTableTag != "LOOKUP_TABLE")
905  {
906  FatalIOErrorInFunction(inFile)
907  << "Expected tag LOOKUP_TABLE but read "
908  << lookupTableTag << nl
909  << exit(FatalIOError);
910  }
911 
912  word lookupTableName(inFile);
913 
914  readField
915  (
916  inFile,
917  selectRegistry(readMode),
918  dataName,
919  dataType,
920  wantedSize//*numComp
921  );
922  }
923  else if (tag == "VECTORS" || tag == "NORMALS")
924  {
925  // 'NORMALS Normals float'
926  string line;
927  inFile.getLine(line);
928  ISpanStream is(line);
929 
930  word dataName(is);
931  word dataType(is);
932  DebugInfo
933  << "Reading vector " << dataName
934  << " of type " << dataType << nl;
935 
936  objectRegistry& reg = selectRegistry(readMode);
937 
938  readField
939  (
940  inFile,
941  reg,
942  dataName,
943  dataType,
944  3*wantedSize
945  );
946 
947  if
948  (
949  vtkDataTypeNames[dataType] == VTK_FLOAT
950  || vtkDataTypeNames[dataType] == VTK_DOUBLE
951  )
952  {
953  objectRegistry::iterator iter = reg.find(dataName);
954  scalarField s(*dynamic_cast<const scalarField*>(iter()));
955  reg.erase(iter);
956  auto fieldVals = autoPtr<vectorIOField>::New
957  (
958  IOobject(dataName, "", reg),
959  s.size()/3
960  );
961 
962  label elemI = 0;
963  for (vector& val : fieldVals())
964  {
965  val.x() = s[elemI++];
966  val.y() = s[elemI++];
967  val.z() = s[elemI++];
968  }
969  regIOobject::store(fieldVals);
970  }
971  }
972  else if (tag == "TEXTURE_COORDINATES")
973  {
974  // 'TEXTURE_COORDINATES TCoords 2 float'
975  string line;
976  inFile.getLine(line);
977  ISpanStream is(line);
978 
979  word dataName(is); //"Tcoords"
980  label dim(readLabel(is));
981  word dataType(is);
982 
983  DebugInfo
984  << "Reading texture coords " << dataName
985  << " dimension " << dim
986  << " of type " << dataType << nl;
987 
988  scalarField coords(dim*points_.size());
989  readBlock(inFile, coords.size(), coords);
990  }
991  else if (tag == "TRIANGLE_STRIPS")
992  {
993  label nStrips(readLabel(inFile));
994  const label nNumbers(readLabel(inFile));
995 
996  labelList elemOffsets, elemVerts;
997 
998  // Total number of faces in strips
999  label nFaces = 0;
1000 
1001  if (version_ < 5.0)
1002  {
1003  // VTK 4.2 and earlier (single-block)
1004  DebugInfo
1005  << "Reading " << nStrips
1006  << " strips (single block)" << nl;
1007 
1008  readBlock(inFile, nNumbers, elemVerts);
1009 
1010  // Count number of faces (triangles)
1011  {
1012  label dataIndex = 0;
1013  for (label i = 0; i < nStrips; ++i)
1014  {
1015  const label nVerts = elemVerts[dataIndex++];
1016  nFaces += nVerts-2;
1017  dataIndex += nVerts;
1018  }
1019  }
1020  }
1021  else
1022  {
1023  // VTK 5.0 and later (OFFSETS and CONNECTIVITY)
1024 
1025  const label nOffsets(nStrips);
1026  --nStrips;
1027 
1028  DebugInfo
1029  << "Reading offsets/connectivity for "
1030  << nStrips << " triangle strips." << nl;
1031 
1032  readOffsetsConnectivity
1033  (
1034  inFile,
1035  "TRIANGLE_STRIPS",
1036  nOffsets, elemOffsets,
1037  nNumbers, elemVerts
1038  );
1039 
1040  // Count number of faces (triangles)
1041  for (label i = 0; i < nStrips; ++i)
1042  {
1043  const label nVerts = (elemOffsets[i+1] - elemOffsets[i]);
1044  nFaces += nVerts-2;
1045  }
1046  }
1047 
1048  // Append into faces
1049  label facei = faces_.size();
1050  faces_.resize(facei+nFaces);
1051  faceMap_.resize(faces_.size());
1052 
1053  label dataIndex = 0;
1054  for (label i = 0; i < nStrips; ++i)
1055  {
1056  const label nVerts =
1057  (
1058  elemOffsets.empty()
1059  ? elemVerts[dataIndex++]
1060  : (elemOffsets[i+1] - elemOffsets[i])
1061  );
1062  const label nTris = nVerts-2;
1063 
1064  // Advance before the first triangle
1065  if (nTris > 0)
1066  {
1067  dataIndex += 2;
1068  }
1069 
1070  for (label triI = 0; triI < nTris; ++triI)
1071  {
1072  faceMap_[facei] = facei;
1073  face& f = faces_[facei++];
1074 
1075  f.resize(3);
1076 
1077  // NOTE: not clear if the orientation is correct
1078  if ((triI % 2) == 0)
1079  {
1080  // Even (eg, 0-1-2, 2-3-4)
1081  f[0] = elemVerts[dataIndex-2];
1082  f[1] = elemVerts[dataIndex-1];
1083  }
1084  else
1085  {
1086  // Odd (eg, 2-1-3, 4-3-5)
1087  f[0] = elemVerts[dataIndex-1];
1088  f[1] = elemVerts[dataIndex-2];
1089  }
1090  f[2] = elemVerts[dataIndex++];
1091  }
1092  }
1093  }
1094  else if (tag == "METADATA")
1095  {
1096  word infoTag(inFile);
1097  if (infoTag != "INFORMATION")
1098  {
1099  FatalIOErrorInFunction(inFile)
1100  << "Unsupported tag "
1101  << infoTag << nl
1102  << exit(FatalIOError);
1103  }
1104  label nInfo(readLabel(inFile));
1105  DebugInfo
1106  << "Ignoring " << nInfo << " metadata information." << nl;
1107 
1108  // Consume rest of line
1109  inFile.getLine(nullptr);
1110  for (label i = 0; i < 2*nInfo; i++)
1111  {
1112  inFile.getLine(nullptr);
1113  }
1114  }
1115  else
1116  {
1117  FatalIOErrorInFunction(inFile)
1118  << "Unsupported tag "
1119  << tag << nl
1120  << exit(FatalIOError);
1121  }
1122  }
1123 
1124 
1125  // There is some problem with orientation of prisms - the point
1126  // ordering seems to be different for some exports (e.g. of cgns)
1127  {
1128  const cellModel& prism = cellModel::ref(cellModel::PRISM);
1129 
1130  label nSwapped = 0;
1131 
1132  for (cellShape& shape : cells_)
1133  {
1134  if (shape.model() == prism)
1135  {
1136  const triPointRef bottom
1137  (
1138  points_[shape[0]],
1139  points_[shape[1]],
1140  points_[shape[2]]
1141  );
1142  const triPointRef top
1143  (
1144  points_[shape[3]],
1145  points_[shape[4]],
1146  points_[shape[5]]
1147  );
1148 
1149  const point bottomCc(bottom.centre());
1150  const vector bottomNormal(bottom.areaNormal());
1151  const point topCc(top.centre());
1152 
1153  if (((topCc - bottomCc) & bottomNormal) < 0)
1154  {
1155  // Flip top and bottom
1156  std::swap(shape[0], shape[3]);
1157  std::swap(shape[1], shape[4]);
1158  std::swap(shape[2], shape[5]);
1159  ++nSwapped;
1160  }
1161  }
1162  }
1163  if (nSwapped)
1164  {
1165  WarningInFunction << "Swapped " << nSwapped
1166  << " prismatic cells" << nl;
1167  }
1168  }
1169 
1170  if (debug)
1171  {
1172  Info<< "Read points:" << points_.size()
1173  << " cells:" << cells_.size()
1174  << " faces:" << faces_.size()
1175  << " lines:" << lines_.size()
1176  << nl << nl;
1177 
1178  Info<< "Cell fields:" << nl;
1179  printFieldStats<vectorIOField>(cellData_);
1180  printFieldStats<scalarIOField>(cellData_);
1181  printFieldStats<labelIOField>(cellData_);
1182  printFieldStats<stringIOList>(cellData_);
1183  Info<< nl << nl;
1184 
1185  Info<< "Point fields:" << nl;
1186  printFieldStats<vectorIOField>(pointData_);
1187  printFieldStats<scalarIOField>(pointData_);
1188  printFieldStats<labelIOField>(pointData_);
1189  printFieldStats<stringIOList>(pointData_);
1190  Info<< nl << nl;
1191 
1192  Info<< "Other fields:" << nl;
1193  printFieldStats<vectorIOField>(otherData_);
1194  printFieldStats<scalarIOField>(otherData_);
1195  printFieldStats<labelIOField>(otherData_);
1196  printFieldStats<stringIOList>(otherData_);
1197  }
1198 }
1199 
1200 
1201 // ************************************************************************* //
SubList< label > labelSubList
A SubList of labels.
Definition: SubList.H:51
static void readBlock(Istream &is, const label n, List< T > &list)
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
type
Types of root.
Definition: Roots.H:52
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::SubStrings splitSpace(const std::string &str, std::string::size_type pos=0)
Split string into sub-strings at whitespace (TAB, NL, VT, FF, CR, SPC)
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:150
IOstream & hex(IOstream &io)
Definition: IOstream.H:560
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:697
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition: label.H:63
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
vtkDataSetType
Enumeration defining the vtk dataset types.
UList< label > labelUList
A UList of labels.
Definition: UList.H:76
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
static const cellModel & ref(const modelType model)
Look up reference to cellModel by enumeration. Fatal on failure.
Definition: cellModels.C:150
static const Enum< parseMode > parseModeNames
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:801
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const labelUList & cellTypes
Definition: setCellMask.H:27
vtkUnstructuredReader(const objectRegistry &obr, ISstream &is)
Construct from input stream, read all.
Vector< scalar > vector
Definition: vector.H:57
#define DebugInfo
Report an information message using Foam::Info.
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
vtkDataType
Enumeration defining the vtk data types.
static bool split(const std::string &line, std::string &key, std::string &val)
Definition: cpuInfo.C:32
List< word > wordList
List of word.
Definition: fileName.H:59
Generic input stream using a standard (STL) stream.
Definition: ISstream.H:51
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:629
static const Enum< vtkDataType > vtkDataTypeNames
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangleFwd.H:39
parseMode
Enumeration defining the parse mode - type of data being read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
static const Enum< vtkDataSetType > vtkDataSetTypeNames
List< label > labelList
A List of labels.
Definition: List.H:61
volScalarField & p
Registry of regIOobjects.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
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))
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...