VTK  9.0.3
vtkEnSightReader.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkEnSightReader.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
20 #ifndef vtkEnSightReader_h
21 #define vtkEnSightReader_h
22 
24 #include "vtkIOEnSightModule.h" // For export macro
25 
26 class vtkDataSet;
28 class vtkEnSightReaderCellIdsType;
29 class vtkIdList;
31 
32 class VTKIOENSIGHT_EXPORT vtkEnSightReader : public vtkGenericEnSightReader
33 {
34 public:
36  void PrintSelf(ostream& os, vtkIndent indent) override;
37 
39  {
40  POINT = 0,
41  BAR2 = 1,
42  BAR3 = 2,
43  NSIDED = 3,
44  TRIA3 = 4,
45  TRIA6 = 5,
46  QUAD4 = 6,
47  QUAD8 = 7,
48  NFACED = 8,
49  TETRA4 = 9,
50  TETRA10 = 10,
51  PYRAMID5 = 11,
52  PYRAMID13 = 12,
53  HEXA8 = 13,
54  HEXA20 = 14,
55  PENTA6 = 15,
56  PENTA15 = 16,
57  NUMBER_OF_ELEMENT_TYPES = 17
58  };
59 
61  {
62  SCALAR_PER_NODE = 0,
63  VECTOR_PER_NODE = 1,
64  TENSOR_SYMM_PER_NODE = 2,
65  SCALAR_PER_ELEMENT = 3,
66  VECTOR_PER_ELEMENT = 4,
67  TENSOR_SYMM_PER_ELEMENT = 5,
68  SCALAR_PER_MEASURED_NODE = 6,
69  VECTOR_PER_MEASURED_NODE = 7,
70  COMPLEX_SCALAR_PER_NODE = 8,
71  COMPLEX_VECTOR_PER_NODE = 9,
72  COMPLEX_SCALAR_PER_ELEMENT = 10,
73  COMPLEX_VECTOR_PER_ELEMENT = 11
74  };
75 
77  {
78  COORDINATES = 0,
79  BLOCK = 1,
80  ELEMENT = 2
81  };
82 
84 
88  vtkGetStringMacro(MeasuredFileName);
90 
92 
96  vtkGetStringMacro(MatchFileName);
98 
99 protected:
101  ~vtkEnSightReader() override;
102 
105 
106  void ClearForNewCaseFileName() override;
107 
109 
112  vtkSetStringMacro(MeasuredFileName);
114 
116 
119  vtkSetStringMacro(MatchFileName);
121 
123 
127  int ReadCaseFileGeometry(char* line);
128  int ReadCaseFileVariable(char* line);
129  int ReadCaseFileTime(char* line);
130  int ReadCaseFileFile(char* line);
132 
133  // set in UpdateInformation to value returned from ReadCaseFile
135 
139  virtual int ReadGeometryFile(
140  const char* fileName, int timeStep, vtkMultiBlockDataSet* output) = 0;
141 
147  const char* fileName, int timeStep, vtkMultiBlockDataSet* output) = 0;
148 
153 
158  virtual int ReadScalarsPerNode(const char* fileName, const char* description, int timeStep,
159  vtkMultiBlockDataSet* output, int measured = 0, int numberOfComponents = 1,
160  int component = 0) = 0;
161 
166  virtual int ReadVectorsPerNode(const char* fileName, const char* description, int timeStep,
167  vtkMultiBlockDataSet* output, int measured = 0) = 0;
168 
173  virtual int ReadTensorsPerNode(
174  const char* fileName, const char* description, int timeStep, vtkMultiBlockDataSet* output) = 0;
175 
180  virtual int ReadScalarsPerElement(const char* fileName, const char* description, int timeStep,
181  vtkMultiBlockDataSet* output, int numberOfComponents = 1, int component = 0) = 0;
182 
188  const char* fileName, const char* description, int timeStep, vtkMultiBlockDataSet* output) = 0;
189 
195  const char* fileName, const char* description, int timeStep, vtkMultiBlockDataSet* output) = 0;
196 
202  int partId, char line[80], const char* name, vtkMultiBlockDataSet* output) = 0;
203 
209  int partId, char line[80], const char* name, vtkMultiBlockDataSet* output) = 0;
210 
214  void AddVariableFileName(const char* fileName1, const char* fileName2 = nullptr);
215 
220 
225 
230  int GetElementType(const char* line);
231 
236  int GetSectionType(const char* line);
237 
241  void ReplaceWildcards(char* filename, int num);
242 
246  void RemoveLeadingBlanks(char* line);
247 
248  // Get the vtkIdList for the given output index and cell type.
249  vtkIdList* GetCellIds(int index, int cellType);
250 
255  void AddToBlock(vtkMultiBlockDataSet* output, unsigned int blockNo, vtkDataSet* dataset);
256 
261  vtkDataSet* GetDataSetFromBlock(vtkMultiBlockDataSet* output, unsigned int blockNo);
262 
266  void SetBlockName(vtkMultiBlockDataSet* output, unsigned int blockNo, const char* name);
267 
269  char* MatchFileName; // may not actually be necessary to read this file
270 
271  // pointer to lists of vtkIdLists (cell ids per element type per part)
272  vtkEnSightReaderCellIdsType* CellIds;
273 
274  // part ids of unstructured outputs
276 
278 
279  // pointers to lists of filenames
280  char** VariableFileNames; // non-complex
282 
283  // array of time sets
286 
287  // array of file sets
290 
291  // collection of filename numbers per time set
294 
295  // collection of filename numbers per file set
298 
299  // collection of number of steps per file per file set
301 
302  // ids of the time and file sets
305 
310 
313 
315  vtkSetMacro(UseTimeSets, vtkTypeBool);
316  vtkGetMacro(UseTimeSets, vtkTypeBool);
317  vtkBooleanMacro(UseTimeSets, vtkTypeBool);
318 
320  vtkSetMacro(UseFileSets, vtkTypeBool);
321  vtkGetMacro(UseFileSets, vtkTypeBool);
322  vtkBooleanMacro(UseFileSets, vtkTypeBool);
323 
325 
326  // global list of points for measured geometry
328 
331 
333 
335 
336 private:
337  vtkEnSightReader(const vtkEnSightReader&) = delete;
338  void operator=(const vtkEnSightReader&) = delete;
339 };
340 
341 #endif
maintain an unordered list of dataset objects
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
superclass for EnSight file readers
int GetSectionType(const char *line)
Determine the section type from a line read a file.
vtkIdList * ComplexVariableFileSetIds
virtual int ReadTensorsPerElement(const char *fileName, const char *description, int timeStep, vtkMultiBlockDataSet *output)=0
Read tensors per element for this dataset.
vtkIdList * UnstructuredPartIds
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkTypeBool UseTimeSets
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
virtual int ReadVectorsPerElement(const char *fileName, const char *description, int timeStep, vtkMultiBlockDataSet *output)=0
Read vectors per element for this dataset.
vtkIdList * FileSets
vtkIdList * GetCellIds(int index, int cellType)
vtkIdListCollection * FileSetFileNameNumbers
void ReplaceWildcards(char *filename, int num)
Replace the *'s in the filename with the given filename number.
vtkTypeBool UseFileSets
void AddToBlock(vtkMultiBlockDataSet *output, unsigned int blockNo, vtkDataSet *dataset)
Convenience method use to convert the readers from VTK 5 multiblock API to the current composite data...
int GetElementType(const char *line)
Determine the element type from a line read a file.
virtual int CreateUnstructuredGridOutput(int partId, char line[80], const char *name, vtkMultiBlockDataSet *output)=0
Read an unstructured part (partId) from the geometry file and create a vtkUnstructuredGrid output.
virtual int ReadVectorsPerNode(const char *fileName, const char *description, int timeStep, vtkMultiBlockDataSet *output, int measured=0)=0
Read vectors per node for this dataset.
void SetBlockName(vtkMultiBlockDataSet *output, unsigned int blockNo, const char *name)
Set the name of a block.
void AddVariableType()
Record the variable type for the variable line just read.
vtkIdList * VariableTimeSetIds
vtkIdListCollection * FileSetNumberOfSteps
void AddVariableDescription(const char *description)
Add another description to the list for a particular variable type.
int ReadCaseFileTime(char *line)
vtkIdList * FileSetsWithFilenameNumbers
vtkEnSightReaderCellIdsType * CellIds
virtual int ReadScalarsPerElement(const char *fileName, const char *description, int timeStep, vtkMultiBlockDataSet *output, int numberOfComponents=1, int component=0)=0
Read scalars per element for this dataset.
int ReadCaseFileVariable(char *line)
int CheckOutputConsistency()
virtual int CreateStructuredGridOutput(int partId, char line[80], const char *name, vtkMultiBlockDataSet *output)=0
Read a structured part from the geometry file and create a vtkStructuredGridOutput.
vtkIdList * VariableFileSetIds
virtual int ReadMeasuredGeometryFile(const char *fileName, int timeStep, vtkMultiBlockDataSet *output)=0
Read the measured geometry file.
int ReadCaseFile()
Read the case file.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
char ** ComplexVariableFileNames
vtkIdList * ComplexVariableTimeSetIds
~vtkEnSightReader() override
vtkIdList * TimeSetsWithFilenameNumbers
virtual int ReadScalarsPerNode(const char *fileName, const char *description, int timeStep, vtkMultiBlockDataSet *output, int measured=0, int numberOfComponents=1, int component=0)=0
Read scalars per node for this dataset.
vtkIdListCollection * TimeSetFileNameNumbers
int ReadCaseFileFile(char *line)
vtkDataSet * GetDataSetFromBlock(vtkMultiBlockDataSet *output, unsigned int blockNo)
Convenience method use to convert the readers from VTK 5 multiblock API to the current composite data...
void ClearForNewCaseFileName() override
Clear data structures such that setting a new case file name works.
int ReadVariableFiles(vtkMultiBlockDataSet *output)
Read the variable files.
int ReadCaseFileGeometry(char *line)
vtkIdList * TimeSetIds
void RemoveLeadingBlanks(char *line)
Remove leading blank spaces from a string.
virtual int ReadGeometryFile(const char *fileName, int timeStep, vtkMultiBlockDataSet *output)=0
Read the geometry file.
virtual int ReadTensorsPerNode(const char *fileName, const char *description, int timeStep, vtkMultiBlockDataSet *output)=0
Read tensors per node for this dataset.
void AddVariableFileName(const char *fileName1, const char *fileName2=nullptr)
Add another file name to the list for a particular variable type.
class to read any type of EnSight files
maintain an ordered list of IdList objects
list of point or cell ids
Definition: vtkIdList.h:31
a simple class to control print indentation
Definition: vtkIndent.h:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Composite dataset that organizes datasets into blocks.
@ component
Definition: vtkX3D.h:181
@ description
Definition: vtkX3D.h:328
@ name
Definition: vtkX3D.h:225
@ index
Definition: vtkX3D.h:252
int vtkTypeBool
Definition: vtkABI.h:69