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