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-2015 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 
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 // Split cellTypes into cells, faces and lines
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[1] = cellVertData[dataIndex++];
306  prismPoints[2] = cellVertData[dataIndex++];
307  prismPoints[3] = cellVertData[dataIndex++];
308  prismPoints[4] = cellVertData[dataIndex++];
309  prismPoints[5] = 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 // Read single field and stores it on the objectRegistry.
361 (
362  ISstream& inFile,
363  objectRegistry& obj,
364  const word& arrayName,
365  const word& dataType,
366  const label size
367 ) const
368 {
369  switch (vtkDataTypeNames[dataType])
370  {
371  case VTK_INT:
372  case VTK_UINT:
373  case VTK_LONG:
374  case VTK_ULONG:
375  case VTK_ID:
376  {
377  autoPtr<labelIOField> fieldVals
378  (
379  new labelIOField
380  (
381  IOobject
382  (
383  arrayName,
384  "",
385  obj
386  ),
387  size
388  )
389  );
390  readBlock(inFile, fieldVals().size(), fieldVals());
391  regIOobject::store(fieldVals);
392  }
393  break;
394 
395  case VTK_FLOAT:
396  case VTK_DOUBLE:
397  {
398  autoPtr<scalarIOField> fieldVals
399  (
400  new scalarIOField
401  (
402  IOobject
403  (
404  arrayName,
405  "",
406  obj
407  ),
408  size
409  )
410  );
411  readBlock(inFile, fieldVals().size(), fieldVals());
412  regIOobject::store(fieldVals);
413  }
414  break;
415 
416  case VTK_STRING:
417  {
418  if (debug)
419  {
420  Info<< "Reading strings:" << size << endl;
421  }
422  autoPtr<stringIOList> fieldVals
423  (
424  new stringIOList
425  (
426  IOobject
427  (
428  arrayName,
429  "",
430  obj
431  ),
432  size
433  )
434  );
435  // Consume current line.
436  inFile.getLine(fieldVals()[0]);
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 // Reads fields, stores them on the objectRegistry. Returns a list of
461 // read fields
463 (
464  ISstream& inFile,
465  objectRegistry& obj,
466  const label wantedSize
467 ) const
468 {
470 
471  word dataName(inFile);
472  if (debug)
473  {
474  Info<< "dataName:" << dataName << endl;
475  }
476  label numArrays(readLabel(inFile));
477  if (debug)
478  {
479  Pout<< "numArrays:" << numArrays << endl;
480  }
481  for (label i = 0; i < numArrays; i++)
482  {
483  word arrayName(inFile);
484  label numComp(readLabel(inFile));
485  label numTuples(readLabel(inFile));
486  word dataType(inFile);
487 
488  if (debug)
489  {
490  Info<< "Reading field " << arrayName
491  << " of " << numTuples << " tuples of rank " << numComp << endl;
492  }
493 
494  if (wantedSize != -1 && numTuples != wantedSize)
495  {
496  FatalIOErrorInFunction(inFile)
497  << "Expected " << wantedSize << " tuples but only have "
498  << numTuples << exit(FatalIOError);
499  }
500 
501  readField
502  (
503  inFile,
504  obj,
505  arrayName,
506  dataType,
507  numTuples*numComp
508  );
509  fields.append(arrayName);
510  }
511  return fields.shrink();
512 }
513 
514 
516 (
517  const parseMode readMode
518 )
519 {
520  if (readMode == CELL_DATA)
521  {
522  return cellData_;
523  }
524  else if (readMode == POINT_DATA)
525  {
526  return pointData_;
527  }
528  else
529  {
530  return otherData_;
531  }
532 }
533 
534 
535 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
536 
538 (
539  const objectRegistry& obr,
540  ISstream& inFile
541 )
542 :
543  cellData_(IOobject("cellData", obr)),
544  pointData_(IOobject("pointData", obr)),
545  otherData_(IOobject("otherData", obr))
546 {
547  read(inFile);
548 }
549 
550 
552 {
553  inFile.getLine(header_);
554  if (debug)
555  {
556  Info<< "Header : " << header_ << endl;
557  }
558  inFile.getLine(title_);
559  if (debug)
560  {
561  Info<< "Title : " << title_ << endl;
562  }
563  inFile.getLine(dataType_);
564  if (debug)
565  {
566  Info<< "dataType : " << dataType_ << endl;
567  }
568 
569  if (dataType_ == "BINARY")
570  {
571  FatalIOErrorInFunction(inFile)
572  << "Binary reading not supported " << exit(FatalIOError);
573  }
574 
575  parseMode readMode = NOMODE;
576  label wantedSize = -1;
577 
578 
579  // Temporary storage for vertices of cells.
580  labelList cellVerts;
581 
582  while (inFile.good())
583  {
584  word tag(inFile);
585 
586  if (!inFile.good())
587  {
588  break;
589  }
590 
591  if (debug)
592  {
593  Info<< "line:" << inFile.lineNumber()
594  << " tag:" << tag << endl;
595  }
596 
597  if (tag == "DATASET")
598  {
599  word geomType(inFile);
600  if (debug)
601  {
602  Info<< "geomType : " << geomType << endl;
603  }
604  readMode = parseModeNames[geomType];
605  wantedSize = -1;
606  }
607  else if (tag == "POINTS")
608  {
609  label nPoints(readLabel(inFile));
610  points_.setSize(nPoints);
611  if (debug)
612  {
613  Info<< "Reading " << nPoints << " numbers representing "
614  << points_.size() << " coordinates." << endl;
615  }
616 
617  word primitiveTag(inFile);
618  if (primitiveTag != "float" && primitiveTag != "double")
619  {
620  FatalIOErrorInFunction(inFile)
621  << "Expected 'float' entry but found "
622  << primitiveTag
623  << exit(FatalIOError);
624  }
625  forAll(points_, i)
626  {
627  inFile >> points_[i].x() >> points_[i].y() >> points_[i].z();
628  }
629  }
630  else if (tag == "CELLS")
631  {
632  label nCells(readLabel(inFile));
633  label nNumbers(readLabel(inFile));
634  if (debug)
635  {
636  Info<< "Reading " << nCells << " cells or faces." << endl;
637  }
638  readBlock(inFile, nNumbers, cellVerts);
639  }
640  else if (tag == "CELL_TYPES")
641  {
642  label nCellTypes(readLabel(inFile));
643 
644  labelList cellTypes;
645  readBlock(inFile, nCellTypes, cellTypes);
646 
647  if (cellTypes.size() > 0 && cellVerts.size() == 0)
648  {
649  FatalIOErrorInFunction(inFile)
650  << "Found " << cellTypes.size()
651  << " cellTypes but no cells."
652  << exit(FatalIOError);
653  }
654 
655  extractCells(inFile, cellTypes, cellVerts);
656  cellVerts.clear();
657  }
658  else if (tag == "LINES")
659  {
660  label nLines(readLabel(inFile));
661  label nNumbers(readLabel(inFile));
662  if (debug)
663  {
664  Info<< "Reading " << nLines << " lines." << endl;
665  }
666  labelList lineVerts;
667  readBlock(inFile, nNumbers, lineVerts);
668 
669  label lineI = lines_.size();
670  lines_.setSize(lineI+nLines);
672 
673  label elemI = 0;
674  for (label i = 0; i < nLines; i++)
675  {
676  lineMap_[lineI] = lineI;
677  labelList& f = lines_[lineI];
678  f.setSize(lineVerts[elemI++]);
679  forAll(f, fp)
680  {
681  f[fp] = lineVerts[elemI++];
682  }
683  lineI++;
684  }
685  }
686  else if (tag == "POLYGONS")
687  {
688  // If in polydata mode
689 
690  label nFaces(readLabel(inFile));
691  label nNumbers(readLabel(inFile));
692  if (debug)
693  {
694  Info<< "Reading " << nFaces << " faces." << endl;
695  }
696  labelList faceVerts;
697  readBlock(inFile, nNumbers, faceVerts);
698 
699  label faceI = faces_.size();
700  faces_.setSize(faceI+nFaces);
702 
703  label elemI = 0;
704  for (label i = 0; i < nFaces; i++)
705  {
706  faceMap_[faceI] = faceI;
707  face& f = faces_[faceI];
708  f.setSize(faceVerts[elemI++]);
709  forAll(f, fp)
710  {
711  f[fp] = faceVerts[elemI++];
712  }
713  faceI++;
714  }
715  }
716  else if (tag == "POINT_DATA")
717  {
718  // 'POINT_DATA 24'
719  readMode = POINT_DATA;
720  wantedSize = points_.size();
721 
722  label nPoints(readLabel(inFile));
723  if (nPoints != wantedSize)
724  {
725  FatalIOErrorInFunction(inFile)
726  << "Reading POINT_DATA : expected " << wantedSize
727  << " but read " << nPoints << exit(FatalIOError);
728  }
729  }
730  else if (tag == "CELL_DATA")
731  {
732  readMode = CELL_DATA;
733  wantedSize = cells_.size()+faces_.size()+lines_.size();
734 
735  label nCells(readLabel(inFile));
736  if (nCells != wantedSize)
737  {
738  FatalIOErrorInFunction(inFile)
739  << "Reading CELL_DATA : expected "
740  << wantedSize
741  << " but read " << nCells << exit(FatalIOError);
742  }
743  }
744  else if (tag == "FIELD")
745  {
746  // wantedSize already set according to type we expected to read.
747  readFieldArray(inFile, selectRegistry(readMode), wantedSize);
748  }
749  else if (tag == "SCALARS")
750  {
751  string line;
752  inFile.getLine(line);
753  IStringStream is(line);
754  word dataName(is);
755  word dataType(is);
756  //label numComp(readLabel(inFile));
757 
758  if (debug)
759  {
760  Info<< "Reading scalar " << dataName
761  << " of type " << dataType
762  << " from lookup table" << endl;
763  }
764 
765  word lookupTableTag(inFile);
766  if (lookupTableTag != "LOOKUP_TABLE")
767  {
768  FatalIOErrorInFunction(inFile)
769  << "Expected tag LOOKUP_TABLE but read "
770  << lookupTableTag
771  << exit(FatalIOError);
772  }
773 
774  word lookupTableName(inFile);
775 
776  readField
777  (
778  inFile,
779  selectRegistry(readMode),
780  dataName,
781  dataType,
782  wantedSize//*numComp
783  );
784  }
785  else if (tag == "VECTORS" || tag == "NORMALS")
786  {
787  // 'NORMALS Normals float'
788  string line;
789  inFile.getLine(line);
790  IStringStream is(line);
791  word dataName(is);
792  word dataType(is);
793  if (debug)
794  {
795  Info<< "Reading vector " << dataName
796  << " of type " << dataType << endl;
797  }
798 
799  objectRegistry& reg = selectRegistry(readMode);
800 
801  readField
802  (
803  inFile,
804  reg,
805  dataName,
806  dataType,
807  3*wantedSize
808  );
809 
810  if
811  (
814  )
815  {
816  objectRegistry::iterator iter = reg.find(dataName);
817  scalarField s(*dynamic_cast<const scalarField*>(iter()));
818  reg.erase(iter);
819  autoPtr<vectorIOField> fieldVals
820  (
821  new vectorIOField
822  (
823  IOobject
824  (
825  dataName,
826  "",
827  reg
828  ),
829  s.size()/3
830  )
831  );
832 
833  label elemI = 0;
834  forAll(fieldVals(), i)
835  {
836  fieldVals()[i].x() = s[elemI++];
837  fieldVals()[i].y() = s[elemI++];
838  fieldVals()[i].z() = s[elemI++];
839  }
840  regIOobject::store(fieldVals);
841  }
842  }
843  else if (tag == "TEXTURE_COORDINATES")
844  {
845  // 'TEXTURE_COORDINATES TCoords 2 float'
846  string line;
847  inFile.getLine(line);
848  IStringStream is(line);
849  word dataName(is); //"Tcoords"
850  label dim(readLabel(is));
851  word dataType(is);
852 
853  if (debug)
854  {
855  Info<< "Reading texture coords " << dataName
856  << " dimension " << dim
857  << " of type " << dataType << endl;
858  }
859 
860  scalarField coords(dim*points_.size());
861  readBlock(inFile, coords.size(), coords);
862  }
863  else if (tag == "TRIANGLE_STRIPS")
864  {
865  label nStrips(readLabel(inFile));
866  label nNumbers(readLabel(inFile));
867  if (debug)
868  {
869  Info<< "Reading " << nStrips << " triangle strips." << endl;
870  }
871  labelList faceVerts;
872  readBlock(inFile, nNumbers, faceVerts);
873 
874  // Count number of triangles
875  label elemI = 0;
876  label nTris = 0;
877  for (label i = 0; i < nStrips; i++)
878  {
879  label nVerts = faceVerts[elemI++];
880  nTris += nVerts-2;
881  elemI += nVerts;
882  }
883 
884 
885  // Store
886  label faceI = faces_.size();
887  faces_.setSize(faceI+nTris);
889  elemI = 0;
890  for (label i = 0; i < nStrips; i++)
891  {
892  label nVerts = faceVerts[elemI++];
893  label nTris = nVerts-2;
894 
895  // Read first triangle
896  faceMap_[faceI] = faceI;
897  face& f = faces_[faceI++];
898  f.setSize(3);
899  f[0] = faceVerts[elemI++];
900  f[1] = faceVerts[elemI++];
901  f[2] = faceVerts[elemI++];
902  for (label triI = 1; triI < nTris; triI++)
903  {
904  faceMap_[faceI] = faceI;
905  face& f = faces_[faceI++];
906  f.setSize(3);
907  f[0] = faceVerts[elemI-1];
908  f[1] = faceVerts[elemI-2];
909  f[2] = faceVerts[elemI++];
910  }
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 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::HashTable< regIOobject * >::iterator
friend class iterator
Declare friendship with the iterator.
Definition: HashTable.H:189
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::vtkUnstructuredReader::header_
string header_
Header.
Definition: vtkUnstructuredReader.H:136
readField
void readField(const IOobject &io, const fvMesh &mesh, const label i, PtrList< GeoField > &fields)
Definition: redistributePar.C:557
Foam::DynamicList< word >
vtkUnstructuredReader.H
Foam::vtkUnstructuredReader::parseModeNames
static const NamedEnum< parseMode, 5 > parseModeNames
Definition: vtkUnstructuredReader.H:107
fields
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
Foam::vtkUnstructuredReader::VTK_FLOAT
@ VTK_FLOAT
Definition: vtkUnstructuredReader.H:76
Foam::vtkUnstructuredReader::readFieldArray
wordList readFieldArray(ISstream &inFile, objectRegistry &obj, const label wantedSize) const
Definition: vtkUnstructuredReader.C:463
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::vtkUnstructuredReader::warnUnhandledType
void warnUnhandledType(Istream &inFile, const label type, labelHashSet &warningGiven) const
Definition: vtkUnstructuredReader.C:86
Foam::vtkUnstructuredReader::vtkDataTypeNames
static const NamedEnum< vtkDataType, 8 > vtkDataTypeNames
Definition: vtkUnstructuredReader.H:82
Foam::ISstream
Generic input stream.
Definition: ISstream.H:51
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::vtkUnstructuredReader::lines_
labelListList lines_
1D cells (=edges)
Definition: vtkUnstructuredReader.H:163
scalarIOField.H
Foam::HashSet< label, Hash< label > >
Foam::vtkUnstructuredReader::vtkDataSetTypeNames
static const NamedEnum< vtkDataSetType, 3 > vtkDataSetTypeNames
Definition: vtkUnstructuredReader.H:93
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::vtkUnstructuredReader::vtkUnstructuredReader
vtkUnstructuredReader(const objectRegistry &obr, ISstream &)
Construct from Istream, read all.
Definition: vtkUnstructuredReader.C:538
cellModeller.H
Foam::vtkUnstructuredReader::read
void read(ISstream &inFile)
Definition: vtkUnstructuredReader.C:551
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::vtkUnstructuredReader::otherData_
objectRegistry otherData_
Other fields.
Definition: vtkUnstructuredReader.H:177
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::ISstream::getLine
ISstream & getLine(string &)
Raw, low-level getline into a string function.
Definition: ISstreamI.H:77
Foam::vtkUnstructuredReader::cellData_
objectRegistry cellData_
Cell based fields.
Definition: vtkUnstructuredReader.H:171
Foam::vtkUnstructuredReader::readField
void readField(ISstream &inFile, objectRegistry &obj, const word &arrayName, const word &dataType, const label size) const
Definition: vtkUnstructuredReader.C:361
Foam::HashTable::erase
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Foam::vtkUnstructuredReader::faceMap_
labelList faceMap_
Map from faces back to original ID.
Definition: vtkUnstructuredReader.H:160
Foam::vtkUnstructuredReader::dataType
const string & dataType() const
DataType.
Definition: vtkUnstructuredReader.H:254
Foam::vtkUnstructuredReader::extractCells
void extractCells(Istream &inFile, const labelList &cellTypes, const labelList &cellVertData)
Definition: vtkUnstructuredReader.C:102
Foam::cellShape
An analytical geometric cellShape.
Definition: cellShape.H:69
Foam::IStringStream
Input from memory buffer stream.
Definition: IStringStream.H:49
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
s
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))
Foam::vtkUnstructuredReader::cells_
cellShapeList cells_
3D cells.
Definition: vtkUnstructuredReader.H:151
Foam::HashTable::find
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Foam::vtkUnstructuredReader::points_
pointField points_
Points.
Definition: vtkUnstructuredReader.H:148
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::hex
IOstream & hex(IOstream &io)
Definition: IOstream.H:564
Foam::tetPoints
Tet storage. Null constructable (unfortunately tetrahedron<point, point> is not)
Definition: tetPoints.H:50
Foam::vtkUnstructuredReader::faces_
faceList faces_
2D cells (=faces)
Definition: vtkUnstructuredReader.H:157
vectorIOField.H
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::vtkUnstructuredReader::lineMap_
labelList lineMap_
Definition: vtkUnstructuredReader.H:165
Foam::regIOobject::store
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:34
Foam::vtkUnstructuredReader::title_
string title_
Title.
Definition: vtkUnstructuredReader.H:139
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
f
labelList f(nPoints)
Foam::vtkUnstructuredReader::pointData_
objectRegistry pointData_
Point based fields.
Definition: vtkUnstructuredReader.H:174
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::vtkUnstructuredReader::dataType_
string dataType_
DataType.
Definition: vtkUnstructuredReader.H:142
Foam::IOList
A List of objects of type <T> with automated input and output.
Definition: IOList.H:50
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:379
Foam::vtkUnstructuredReader::readBlock
void readBlock(Istream &inFile, const label n, List< T > &lst) const
Definition: vtkUnstructuredReaderTemplates.C:37
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::line
A line primitive.
Definition: line.H:56
Foam::IOstream::lineNumber
label lineNumber() const
Return current stream line number.
Definition: IOstream.H:438
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::cellModel
Maps a geometry to a set of cell primitives, which enables geometric cell data to be calculated witho...
Definition: cellModel.H:64
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
IOWarningInFunction
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Definition: messageStream.H:271
Foam::vtkUnstructuredReader::NOMODE
@ NOMODE
Definition: vtkUnstructuredReader.H:100
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::vtkUnstructuredReader::POINT_DATA
@ POINT_DATA
Definition: vtkUnstructuredReader.H:104
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::vtkUnstructuredReader::CELL_DATA
@ CELL_DATA
Definition: vtkUnstructuredReader.H:103
Foam::vtkUnstructuredReader::selectRegistry
objectRegistry & selectRegistry(const parseMode readMode)
Definition: vtkUnstructuredReader.C:516
Foam::vtkUnstructuredReader::VTK_DOUBLE
@ VTK_DOUBLE
Definition: vtkUnstructuredReader.H:77
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Foam::vtkUnstructuredReader::parseMode
parseMode
Enumeration defining the parse mode - what type of data is being.
Definition: vtkUnstructuredReader.H:98
labelIOField.H
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
stringIOList.H