vtkUnstructuredReader.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2019 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 \*---------------------------------------------------------------------------*/
25 
26 #include "vtkUnstructuredReader.H"
27 #include "labelIOField.H"
28 #include "scalarIOField.H"
29 #include "stringIOList.H"
30 #include "cellModeller.H"
31 #include "vectorIOField.H"
32 #include "stringOps.H"
33 
34 /* * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * */
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(vtkUnstructuredReader, 1);
39 
40  template<>
41  const char*
43  {
44  "int",
45  "unsigned_int",
46  "long",
47  "unsigned_long",
48  "float",
49  "double",
50  "string",
51  "vtkIdType"
52  };
55 
56 
57  template<>
58  const char*
60  {
61  "FIELD",
62  "SCALARS",
63  "VECTORS"
64  };
67 
68 
69  template<>
70  const char*
72  {
73  "NOMODE",
74  "UNSTRUCTURED_GRID",
75  "POLYDATA",
76  "CELL_DATA",
77  "POINT_DATA"
78  };
81 }
82 
83 
84 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
85 
86 void Foam::vtkUnstructuredReader::warnUnhandledType
87 (
88  Istream& inFile,
89  const label type,
90  labelHashSet& warningGiven
91 ) const
92 {
93  if (warningGiven.insert(type))
94  {
95  IOWarningInFunction(inFile)
96  << "Skipping unknown cell type " << type << endl;
97  }
98 }
99 
100 
101 void Foam::vtkUnstructuredReader::extractCells
102 (
103  Istream& inFile,
104  const labelList& cellTypes,
105  const labelList& cellVertData
106 )
107 {
108  const cellModel& hex = *(cellModeller::lookup("hex"));
109  const cellModel& prism = *(cellModeller::lookup("prism"));
110  const cellModel& pyr = *(cellModeller::lookup("pyr"));
111  const cellModel& tet = *(cellModeller::lookup("tet"));
112 
113  labelList tetPoints(4);
114  labelList pyrPoints(5);
115  labelList prismPoints(6);
116  labelList hexPoints(8);
117 
118  label celli = cells_.size();
119  cells_.setSize(celli+cellTypes.size());
120  cellMap_.setSize(cells_.size(), -1);
121 
122  label facei = faces_.size();
123  faces_.setSize(facei+cellTypes.size());
124  faceMap_.setSize(faces_.size(), -1);
125 
126  label lineI = lines_.size();
127  lines_.setSize(lineI+cellTypes.size());
128  lineMap_.setSize(lines_.size(), -1);
129 
130  label dataIndex = 0;
131 
132 
133  // To mark whether unhandled type has been visited.
134  labelHashSet warningGiven;
135 
136  forAll(cellTypes, i)
137  {
138  switch (cellTypes[i])
139  {
140  case VTK_VERTEX:
141  {
142  warnUnhandledType(inFile, cellTypes[i], warningGiven);
143  label nRead = cellVertData[dataIndex++];
144  if (nRead != 1)
145  {
147  (
148  inFile
149  ) << "Expected size 1 for VTK_VERTEX but found "
150  << nRead << exit(FatalIOError);
151  }
152  dataIndex += nRead;
153  }
154  break;
155 
156  case VTK_POLY_VERTEX:
157  {
158  warnUnhandledType(inFile, cellTypes[i], warningGiven);
159  label nRead = cellVertData[dataIndex++];
160  dataIndex += nRead;
161  }
162  break;
163 
164  case VTK_LINE:
165  {
166  // warnUnhandledType(inFile, cellTypes[i], warningGiven);
167  label nRead = cellVertData[dataIndex++];
168  if (nRead != 2)
169  {
171  (
172  inFile
173  ) << "Expected size 2 for VTK_LINE but found "
174  << nRead << exit(FatalIOError);
175  }
176  lineMap_[lineI] = i;
177  labelList& segment = lines_[lineI++];
178  segment.setSize(2);
179  segment[0] = cellVertData[dataIndex++];
180  segment[1] = cellVertData[dataIndex++];
181  }
182  break;
183 
184  case VTK_POLY_LINE:
185  {
186  // warnUnhandledType(inFile, cellTypes[i], warningGiven);
187  label nRead = cellVertData[dataIndex++];
188  lineMap_[lineI] = i;
189  labelList& segment = lines_[lineI++];
190  segment.setSize(nRead);
191  forAll(segment, i)
192  {
193  segment[i] = cellVertData[dataIndex++];
194  }
195  }
196  break;
197 
198  case VTK_TRIANGLE:
199  {
200  faceMap_[facei] = i;
201  face& f = faces_[facei++];
202  f.setSize(3);
203  label nRead = cellVertData[dataIndex++];
204  if (nRead != 3)
205  {
207  (
208  inFile
209  ) << "Expected size 3 for VTK_TRIANGLE but found "
210  << nRead << exit(FatalIOError);
211  }
212  f[0] = cellVertData[dataIndex++];
213  f[1] = cellVertData[dataIndex++];
214  f[2] = cellVertData[dataIndex++];
215  }
216  break;
217 
218  case VTK_QUAD:
219  {
220  faceMap_[facei] = i;
221  face& f = faces_[facei++];
222  f.setSize(4);
223  label nRead = cellVertData[dataIndex++];
224  if (nRead != 4)
225  {
227  (
228  inFile
229  ) << "Expected size 4 for VTK_QUAD but found "
230  << nRead << exit(FatalIOError);
231  }
232  f[0] = cellVertData[dataIndex++];
233  f[1] = cellVertData[dataIndex++];
234  f[2] = cellVertData[dataIndex++];
235  f[3] = cellVertData[dataIndex++];
236  }
237  break;
238 
239  case VTK_POLYGON:
240  {
241  faceMap_[facei] = i;
242  face& f = faces_[facei++];
243  label nRead = cellVertData[dataIndex++];
244  f.setSize(nRead);
245  forAll(f, fp)
246  {
247  f[fp] = cellVertData[dataIndex++];
248  }
249  }
250  break;
251 
252  case VTK_TETRA:
253  {
254  label nRead = cellVertData[dataIndex++];
255  if (nRead != 4)
256  {
258  (
259  inFile
260  ) << "Expected size 4 for VTK_TETRA but found "
261  << nRead << exit(FatalIOError);
262  }
263  tetPoints[0] = cellVertData[dataIndex++];
264  tetPoints[1] = cellVertData[dataIndex++];
265  tetPoints[2] = cellVertData[dataIndex++];
266  tetPoints[3] = cellVertData[dataIndex++];
267  cellMap_[celli] = i;
268  cells_[celli++] = cellShape(tet, tetPoints, true);
269  }
270  break;
271 
272  case VTK_PYRAMID:
273  {
274  label nRead = cellVertData[dataIndex++];
275  if (nRead != 5)
276  {
278  (
279  inFile
280  ) << "Expected size 5 for VTK_PYRAMID but found "
281  << nRead << exit(FatalIOError);
282  }
283  pyrPoints[0] = cellVertData[dataIndex++];
284  pyrPoints[1] = cellVertData[dataIndex++];
285  pyrPoints[2] = cellVertData[dataIndex++];
286  pyrPoints[3] = cellVertData[dataIndex++];
287  pyrPoints[4] = cellVertData[dataIndex++];
288  cellMap_[celli] = i;
289  cells_[celli++] = cellShape(pyr, pyrPoints, true);
290  }
291  break;
292 
293  case VTK_WEDGE:
294  {
295  label nRead = cellVertData[dataIndex++];
296  if (nRead != 6)
297  {
299  (
300  inFile
301  ) << "Expected size 6 for VTK_WEDGE but found "
302  << nRead << exit(FatalIOError);
303  }
304  prismPoints[0] = cellVertData[dataIndex++];
305  prismPoints[2] = cellVertData[dataIndex++];
306  prismPoints[1] = cellVertData[dataIndex++];
307  prismPoints[3] = cellVertData[dataIndex++];
308  prismPoints[5] = cellVertData[dataIndex++];
309  prismPoints[4] = cellVertData[dataIndex++];
310  cellMap_[celli] = i;
311  cells_[celli++] = cellShape(prism, prismPoints, true);
312  }
313  break;
314 
315  case VTK_HEXAHEDRON:
316  {
317  label nRead = cellVertData[dataIndex++];
318  if (nRead != 8)
319  {
321  (
322  inFile
323  ) << "Expected size 8 for VTK_HEXAHEDRON but found "
324  << nRead << exit(FatalIOError);
325  }
326  hexPoints[0] = cellVertData[dataIndex++];
327  hexPoints[1] = cellVertData[dataIndex++];
328  hexPoints[2] = cellVertData[dataIndex++];
329  hexPoints[3] = cellVertData[dataIndex++];
330  hexPoints[4] = cellVertData[dataIndex++];
331  hexPoints[5] = cellVertData[dataIndex++];
332  hexPoints[6] = cellVertData[dataIndex++];
333  hexPoints[7] = cellVertData[dataIndex++];
334  cellMap_[celli] = i;
335  cells_[celli++] = cellShape(hex, hexPoints, true);
336  }
337  break;
338 
339  default:
340  warnUnhandledType(inFile, cellTypes[i], warningGiven);
341  label nRead = cellVertData[dataIndex++];
342  dataIndex += nRead;
343  }
344  }
345 
346  if (debug)
347  {
348  Info<< "Read " << celli << " cells;" << facei << " faces." << endl;
349  }
350  cells_.setSize(celli);
351  cellMap_.setSize(celli);
352  faces_.setSize(facei);
353  faceMap_.setSize(facei);
354  lines_.setSize(lineI);
355  lineMap_.setSize(lineI);
356 }
357 
358 
359 void Foam::vtkUnstructuredReader::readField
360 (
361  ISstream& inFile,
362  objectRegistry& obj,
363  const word& arrayName,
364  const word& dataType,
365  const label size
366 ) const
367 {
368  switch (vtkDataTypeNames[dataType])
369  {
370  case VTK_INT:
371  case VTK_UINT:
372  case VTK_LONG:
373  case VTK_ULONG:
374  case VTK_ID:
375  {
376  autoPtr<labelIOField> fieldVals
377  (
378  new labelIOField
379  (
380  IOobject
381  (
382  arrayName,
383  "",
384  obj
385  ),
386  size
387  )
388  );
389  readBlock(inFile, fieldVals().size(), fieldVals());
390  regIOobject::store(fieldVals);
391  }
392  break;
393 
394  case VTK_FLOAT:
395  case VTK_DOUBLE:
396  {
397  autoPtr<scalarIOField> fieldVals
398  (
399  new scalarIOField
400  (
401  IOobject
402  (
403  arrayName,
404  "",
405  obj
406  ),
407  size
408  )
409  );
410  readBlock(inFile, fieldVals().size(), fieldVals());
411  regIOobject::store(fieldVals);
412  }
413  break;
414 
415  case VTK_STRING:
416  {
417  if (debug)
418  {
419  Info<< "Reading strings:" << size << endl;
420  }
421  autoPtr<stringIOList> fieldVals
422  (
423  new stringIOList
424  (
425  IOobject
426  (
427  arrayName,
428  "",
429  obj
430  ),
431  size
432  )
433  );
434  // Consume current line.
435  inFile.getLine(fieldVals()[0]);
436 
437  // Read without parsing
438  forAll(fieldVals(), i)
439  {
440  inFile.getLine(fieldVals()[i]);
441  }
442  regIOobject::store(fieldVals);
443  }
444  break;
445 
446  default:
447  {
448  IOWarningInFunction(inFile)
449  << "Unhandled type " << vtkDataTypeNames[dataType] << endl
450  << "Skipping " << size
451  << " words." << endl;
452  scalarField fieldVals;
453  readBlock(inFile, size, fieldVals);
454  }
455  break;
456  }
457 }
458 
459 
460 Foam::wordList Foam::vtkUnstructuredReader::readFieldArray
461 (
462  ISstream& inFile,
463  objectRegistry& obj,
464  const label wantedSize
465 ) const
466 {
468 
469  word dataName(inFile);
470  if (debug)
471  {
472  Info<< "dataName:" << dataName << endl;
473  }
474  label numArrays(readLabel(inFile));
475  if (debug)
476  {
477  Pout<< "numArrays:" << numArrays << endl;
478  }
479  for (label i = 0; i < numArrays; i++)
480  {
481  word arrayName(inFile);
482  label numComp(readLabel(inFile));
483  label numTuples(readLabel(inFile));
484  word dataType(inFile);
485 
486  if (debug)
487  {
488  Info<< "Reading field " << arrayName
489  << " of " << numTuples << " tuples of rank " << numComp << endl;
490  }
491 
492  if (wantedSize != -1 && numTuples != wantedSize)
493  {
494  FatalIOErrorInFunction(inFile)
495  << "Expected " << wantedSize << " tuples but only have "
496  << numTuples << exit(FatalIOError);
497  }
498 
499  readField
500  (
501  inFile,
502  obj,
503  arrayName,
504  dataType,
505  numTuples*numComp
506  );
507  fields.append(arrayName);
508  }
509  return fields.shrink();
510 }
511 
512 
513 Foam::objectRegistry& Foam::vtkUnstructuredReader::selectRegistry
514 (
515  const parseMode readMode
516 )
517 {
518  if (readMode == CELL_DATA)
519  {
520  return cellData_;
521  }
522  else if (readMode == POINT_DATA)
523  {
524  return pointData_;
525  }
526  else
527  {
528  return otherData_;
529  }
530 }
531 
532 
533 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
534 
536 (
537  const objectRegistry& obr,
538  ISstream& inFile
539 )
540 :
541  cellData_(IOobject("cellData", obr)),
542  pointData_(IOobject("pointData", obr)),
543  otherData_(IOobject("otherData", obr))
544 {
545  read(inFile);
546 }
547 
548 
549 void Foam::vtkUnstructuredReader::read(ISstream& inFile)
550 {
551  inFile.getLine(header_);
552  if (debug)
553  {
554  Info<< "Header : " << header_ << endl;
555  }
556  inFile.getLine(title_);
557  if (debug)
558  {
559  Info<< "Title : " << title_ << endl;
560  }
561  inFile.getLine(dataType_);
562  if (debug)
563  {
564  Info<< "dataType : " << dataType_ << endl;
565  }
566 
567  if (dataType_ == "BINARY")
568  {
569  FatalIOErrorInFunction(inFile)
570  << "Binary reading not supported " << exit(FatalIOError);
571  }
572 
573  parseMode readMode = NOMODE;
574  label wantedSize = -1;
575 
576 
577  // Temporary storage for vertices of cells.
578  labelList cellVerts;
579 
580  while (inFile.good())
581  {
582  word tag(inFile);
583 
584  if (!inFile.good())
585  {
586  break;
587  }
588 
589  if (debug)
590  {
591  Info<< "line:" << inFile.lineNumber()
592  << " tag:" << tag << endl;
593  }
594 
595  if (tag == "DATASET")
596  {
597  word geomType(inFile);
598  if (debug)
599  {
600  Info<< "geomType : " << geomType << endl;
601  }
602  readMode = parseModeNames[geomType];
603  wantedSize = -1;
604  }
605  else if (tag == "POINTS")
606  {
607  label nPoints(readLabel(inFile));
608  points_.setSize(nPoints);
609  if (debug)
610  {
611  Info<< "Reading " << nPoints << " numbers representing "
612  << points_.size() << " coordinates." << endl;
613  }
614 
615  word primitiveTag(inFile);
616  if (primitiveTag != "float" && primitiveTag != "double")
617  {
618  FatalIOErrorInFunction(inFile)
619  << "Expected 'float' entry but found "
620  << primitiveTag
621  << exit(FatalIOError);
622  }
623  forAll(points_, i)
624  {
625  inFile >> points_[i].x() >> points_[i].y() >> points_[i].z();
626  }
627  }
628  else if (tag == "CELLS")
629  {
630  label nCells(readLabel(inFile));
631  label nNumbers(readLabel(inFile));
632  if (debug)
633  {
634  Info<< "Reading " << nCells << " cells or faces." << endl;
635  }
636  readBlock(inFile, nNumbers, cellVerts);
637  }
638  else if (tag == "CELL_TYPES")
639  {
640  label nCellTypes(readLabel(inFile));
641 
642  labelList cellTypes;
643  readBlock(inFile, nCellTypes, cellTypes);
644 
645  if (cellTypes.size() > 0 && cellVerts.size() == 0)
646  {
647  FatalIOErrorInFunction(inFile)
648  << "Found " << cellTypes.size()
649  << " cellTypes but no cells."
650  << exit(FatalIOError);
651  }
652 
653  extractCells(inFile, cellTypes, cellVerts);
654  cellVerts.clear();
655  }
656  else if (tag == "LINES")
657  {
658  label nLines(readLabel(inFile));
659  label nNumbers(readLabel(inFile));
660  if (debug)
661  {
662  Info<< "Reading " << nLines << " lines." << endl;
663  }
664  labelList lineVerts;
665  readBlock(inFile, nNumbers, lineVerts);
666 
667  label lineI = lines_.size();
668  lines_.setSize(lineI+nLines);
669  lineMap_.setSize(lines_.size());
670 
671  label elemI = 0;
672  for (label i = 0; i < nLines; i++)
673  {
674  lineMap_[lineI] = lineI;
675  labelList& f = lines_[lineI];
676  f.setSize(lineVerts[elemI++]);
677  forAll(f, fp)
678  {
679  f[fp] = lineVerts[elemI++];
680  }
681  lineI++;
682  }
683  }
684  else if (tag == "POLYGONS")
685  {
686  // If in polydata mode
687 
688  label nFaces(readLabel(inFile));
689  label nNumbers(readLabel(inFile));
690  if (debug)
691  {
692  Info<< "Reading " << nFaces << " faces." << endl;
693  }
694  labelList faceVerts;
695  readBlock(inFile, nNumbers, faceVerts);
696 
697  label facei = faces_.size();
698  faces_.setSize(facei+nFaces);
699  faceMap_.setSize(faces_.size());
700 
701  label elemI = 0;
702  for (label i = 0; i < nFaces; i++)
703  {
704  faceMap_[facei] = facei;
705  face& f = faces_[facei];
706  f.setSize(faceVerts[elemI++]);
707  forAll(f, fp)
708  {
709  f[fp] = faceVerts[elemI++];
710  }
711  facei++;
712  }
713  }
714  else if (tag == "POINT_DATA")
715  {
716  // 'POINT_DATA 24'
717  readMode = POINT_DATA;
718  wantedSize = points_.size();
719 
720  label nPoints(readLabel(inFile));
721  if (nPoints != wantedSize)
722  {
723  FatalIOErrorInFunction(inFile)
724  << "Reading POINT_DATA : expected " << wantedSize
725  << " but read " << nPoints << exit(FatalIOError);
726  }
727  }
728  else if (tag == "CELL_DATA")
729  {
730  readMode = CELL_DATA;
731  wantedSize = cells_.size()+faces_.size()+lines_.size();
732 
733  label nCells(readLabel(inFile));
734  if (nCells != wantedSize)
735  {
736  FatalIOErrorInFunction(inFile)
737  << "Reading CELL_DATA : expected "
738  << wantedSize
739  << " but read " << nCells << exit(FatalIOError);
740  }
741  }
742  else if (tag == "FIELD")
743  {
744  // wantedSize already set according to type we expected to read.
745  readFieldArray(inFile, selectRegistry(readMode), wantedSize);
746  }
747  else if (tag == "SCALARS")
748  {
749  string line;
750  inFile.getLine(line);
751  IStringStream is(line);
752  word dataName(is);
753  word dataType(is);
754  // label numComp(readLabel(inFile));
755 
756  if (debug)
757  {
758  Info<< "Reading scalar " << dataName
759  << " of type " << dataType
760  << " from lookup table" << endl;
761  }
762 
763  word lookupTableTag(inFile);
764  if (lookupTableTag != "LOOKUP_TABLE")
765  {
766  FatalIOErrorInFunction(inFile)
767  << "Expected tag LOOKUP_TABLE but read "
768  << lookupTableTag
769  << exit(FatalIOError);
770  }
771 
772  word lookupTableName(inFile);
773 
774  readField
775  (
776  inFile,
777  selectRegistry(readMode),
778  dataName,
779  dataType,
780  wantedSize//*numComp
781  );
782  }
783  else if (tag == "VECTORS" || tag == "NORMALS")
784  {
785  // 'NORMALS Normals float'
786  string line;
787  inFile.getLine(line);
788  IStringStream is(line);
789  word dataName(is);
790  word dataType(is);
791  if (debug)
792  {
793  Info<< "Reading vector " << dataName
794  << " of type " << dataType << endl;
795  }
796 
797  objectRegistry& reg = selectRegistry(readMode);
798 
799  readField
800  (
801  inFile,
802  reg,
803  dataName,
804  dataType,
805  3*wantedSize
806  );
807 
808  if
809  (
810  vtkDataTypeNames[dataType] == VTK_FLOAT
811  || vtkDataTypeNames[dataType] == VTK_DOUBLE
812  )
813  {
814  objectRegistry::iterator iter = reg.find(dataName);
815  scalarField s(*dynamic_cast<const scalarField*>(iter()));
816  reg.erase(iter);
817  autoPtr<vectorIOField> fieldVals
818  (
819  new vectorIOField
820  (
821  IOobject
822  (
823  dataName,
824  "",
825  reg
826  ),
827  s.size()/3
828  )
829  );
830 
831  label elemI = 0;
832  forAll(fieldVals(), i)
833  {
834  fieldVals()[i].x() = s[elemI++];
835  fieldVals()[i].y() = s[elemI++];
836  fieldVals()[i].z() = s[elemI++];
837  }
838  regIOobject::store(fieldVals);
839  }
840  }
841  else if (tag == "TEXTURE_COORDINATES")
842  {
843  // 'TEXTURE_COORDINATES TCoords 2 float'
844  string line;
845  inFile.getLine(line);
846  IStringStream is(line);
847  word dataName(is); //"Tcoords"
848  label dim(readLabel(is));
849  word dataType(is);
850 
851  if (debug)
852  {
853  Info<< "Reading texture coords " << dataName
854  << " dimension " << dim
855  << " of type " << dataType << endl;
856  }
857 
858  scalarField coords(dim*points_.size());
859  readBlock(inFile, coords.size(), coords);
860  }
861  else if (tag == "TRIANGLE_STRIPS")
862  {
863  label nStrips(readLabel(inFile));
864  label nNumbers(readLabel(inFile));
865  if (debug)
866  {
867  Info<< "Reading " << nStrips << " triangle strips." << endl;
868  }
869  labelList faceVerts;
870  readBlock(inFile, nNumbers, faceVerts);
871 
872  // Count number of triangles
873  label elemI = 0;
874  label nTris = 0;
875  for (label i = 0; i < nStrips; i++)
876  {
877  label nVerts = faceVerts[elemI++];
878  nTris += nVerts-2;
879  elemI += nVerts;
880  }
881 
882 
883  // Store
884  label facei = faces_.size();
885  faces_.setSize(facei+nTris);
886  faceMap_.setSize(faces_.size());
887  elemI = 0;
888  for (label i = 0; i < nStrips; i++)
889  {
890  label nVerts = faceVerts[elemI++];
891  label nTris = nVerts-2;
892 
893  // Read first triangle
894  faceMap_[facei] = facei;
895  face& f = faces_[facei++];
896  f.setSize(3);
897  f[0] = faceVerts[elemI++];
898  f[1] = faceVerts[elemI++];
899  f[2] = faceVerts[elemI++];
900  for (label triI = 1; triI < nTris; triI++)
901  {
902  faceMap_[facei] = facei;
903  face& f = faces_[facei++];
904  f.setSize(3);
905  f[0] = faceVerts[elemI-1];
906  f[1] = faceVerts[elemI-2];
907  f[2] = faceVerts[elemI++];
908  }
909  }
910  }
911  else if (tag == "METADATA")
912  {
913  // Read rest of the line
914  string line;
915  inFile.getLine(line);
916 
917  // Skip until an empty line is found
918  inFile.getLine(line);
919  while (!stringOps::trim(line).empty())
920  {
921  inFile.getLine(line);
922  }
923  }
924  else
925  {
926  FatalIOErrorInFunction(inFile)
927  << "Unsupported tag "
928  << tag << exit(FatalIOError);
929  }
930  }
931 
932  if (debug)
933  {
934  Info<< "Read points:" << points_.size()
935  << " cellShapes:" << cells_.size()
936  << " faces:" << faces_.size()
937  << " lines:" << lines_.size()
938  << nl << endl;
939 
940  Info<< "Cell fields:" << endl;
941  printFieldStats<vectorIOField>(cellData_);
942  printFieldStats<scalarIOField>(cellData_);
943  printFieldStats<labelIOField>(cellData_);
944  printFieldStats<stringIOList>(cellData_);
945  Info<< nl << endl;
946 
947  Info<< "Point fields:" << endl;
948  printFieldStats<vectorIOField>(pointData_);
949  printFieldStats<scalarIOField>(pointData_);
950  printFieldStats<labelIOField>(pointData_);
951  printFieldStats<stringIOList>(pointData_);
952  Info<< nl << endl;
953 
954  Info<< "Other fields:" << endl;
955  printFieldStats<vectorIOField>(otherData_);
956  printFieldStats<scalarIOField>(otherData_);
957  printFieldStats<labelIOField>(otherData_);
958  printFieldStats<stringIOList>(otherData_);
959  }
960 }
961 
962 
963 // ************************************************************************* //
vtkUnstructuredReader(const objectRegistry &obr, ISstream &)
Construct from Istream, read all.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
static const cellModel * lookup(const word &)
Look up a model by name and return a pointer to the model or nullptr.
Definition: cellModeller.C:100
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
IOstream & hex(IOstream &io)
Definition: IOstream.H:564
An analytical geometric cellShape.
Definition: cellShape.H:69
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A List of objects of type <T> with automated input and output.
Definition: IOList.H:50
label lineNumber() const
Return current stream line number.
Definition: IOstream.H:438
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:108
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
static const NamedEnum< vtkDataSetType, 3 > vtkDataSetTypeNames
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:371
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
string trim(const string &)
Return string trimmed of leading and trailing whitespace.
Definition: stringOps.C:923
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
bool read(const char *, int32_t &)
Definition: int32IO.C:85
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))
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:124
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
label readLabel(Istream &is)
Definition: label.H:64
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:34
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
static const char nl
Definition: Ostream.H:265
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
void setSize(const label)
Reset size of List.
Definition: List.C:281
Generic input stream.
Definition: ISstream.H:51
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
Input from memory buffer stream.
Definition: IStringStream.H:49
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Maps a geometry to a set of cell primitives, which enables geometric cell data to be calculated witho...
Definition: cellModel.H:64
parseMode
Enumeration defining the parse mode - what type of data is being.
messageStream Info
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Registry of regIOobjects.
friend class iterator
Declare friendship with the iterator.
Definition: HashTable.H:194
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static const NamedEnum< vtkDataType, 8 > vtkDataTypeNames
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
Namespace for OpenFOAM.
static const NamedEnum< parseMode, 5 > parseModeNames
ISstream & getLine(string &)
Raw, low-level getline into a string function.
Definition: ISstreamI.H:77
IOerror FatalIOError