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