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