VTK  9.3.20240327
vtkTecplotReader.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright (c) 2000 - 2009, Lawrence Livermore National Security, LLC
3 // SPDX-License-Identifier: BSD-3-Clause
4 
54 #ifndef vtkTecplotReader_h
55 #define vtkTecplotReader_h
56 
57 #include "vtkIOGeometryModule.h" // For export macro
59 
60 #include <string> // STL Header; Required for string
61 #include <vector> // STL Header; Required for vector
62 
63 VTK_ABI_NAMESPACE_BEGIN
64 class vtkPoints;
65 class vtkCellData;
66 class vtkPointData;
67 class vtkCallbackCommand;
71 class vtkTecplotReaderInternal;
72 
73 class VTKIOGEOMETRY_EXPORT vtkTecplotReader : public vtkMultiBlockDataSetAlgorithm
74 {
75 public:
76  static vtkTecplotReader* New();
78  void PrintSelf(ostream& os, vtkIndent indent) override;
79 
81 
84  vtkGetMacro(NumberOfVariables, int);
86 
90  void SetFileName(VTK_FILEPATH const char* fileName);
91 
95  const char* GetDataTitle();
96 
101 
106  const char* GetBlockName(int blockIdx);
107 
113 
118  const char* GetDataAttributeName(int attrIndx);
119 
125  int IsDataAttributeCellBased(const char* attrName);
126 
132  int IsDataAttributeCellBased(int attrIndx);
133 
138 
142  const char* GetDataArrayName(int arrayIdx);
143 
147  int GetDataArrayStatus(const char* arayName);
148 
153  void SetDataArrayStatus(const char* arayName, int bChecked);
154 
155 protected:
157  ~vtkTecplotReader() override;
158 
161  vtkInformationVector* outputVector) override;
163 
167  static void SelectionModifiedCallback(vtkObject*, unsigned long, void* tpReader, void*);
168 
174  void Init();
175 
180 
186 
194  int numNodes, int numCells, vtkPoints* theNodes, vtkPointData* nodeData, vtkCellData* cellData);
195 
204  void GetArraysFromPointPackingZone(int numNodes, vtkPoints* theNodes, vtkPointData* nodeData);
205 
213  void GetStructuredGridFromBlockPackingZone(int iDimSize, int jDimSize, int kDimSize, int zoneIndx,
214  const char* zoneName, vtkMultiBlockDataSet* multZone);
215 
223  void GetStructuredGridFromPointPackingZone(int iDimSize, int jDimSize, int kDimSize, int zoneIndx,
224  const char* zoneName, vtkMultiBlockDataSet* multZone);
225 
233  void GetUnstructuredGridFromBlockPackingZone(int numNodes, int numCells, const char* cellType,
234  int zoneIndx, const char* zoneName, vtkMultiBlockDataSet* multZone);
235 
243  void GetPolyhedralGridFromBlockPackingZone(int numNodes, int numElements, int numFaces,
244  int zoneIndex, const char* zoneName, vtkMultiBlockDataSet* multZone);
245 
253  void GetPolygonalGridFromBlockPackingZone(int numNodes, int numElements, int numFaces,
254  int zoneIndex, const char* zoneName, vtkMultiBlockDataSet* multZone);
255 
260  void GetPolyhedralGridCells(int numberCells, int numFaces, vtkUnstructuredGrid* unstruct) const;
261 
266  void GetPolygonalGridCells(int numFaces, int numEdges, vtkUnstructuredGrid* unstruct) const;
267 
275  void GetUnstructuredGridFromPointPackingZone(int numNodes, int numCells, const char* cellType,
276  int zoneIndx, const char* zoneName, vtkMultiBlockDataSet* multZone);
277 
283  int numberCells, const char* cellTypeStr, vtkUnstructuredGrid* unstrctGrid);
284 
286  char* FileName;
289  vtkTecplotReaderInternal* Internal;
290 
292  std::vector<int> CellBased;
293  std::vector<std::string> ZoneNames;
294  std::vector<std::string> Variables;
295 
296 private:
297  vtkTecplotReader(const vtkTecplotReader&) = delete;
298  void operator=(const vtkTecplotReader&) = delete;
299 };
300 
301 VTK_ABI_NAMESPACE_END
302 #endif
supports function callbacks
represent and manipulate cell attribute data
Definition: vtkCellData.h:140
Store on/off settings for data arrays, etc.
a simple class to control print indentation
Definition: vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
Composite dataset that organizes datasets into blocks.
abstract base class for most VTK objects
Definition: vtkObject.h:161
represent and manipulate point attribute data
Definition: vtkPointData.h:139
represent and manipulate 3D points
Definition: vtkPoints.h:138
A concrete class to read an ASCII Tecplot file.
void ReadFile(vtkMultiBlockDataSet *multZone)
This function, the data loading engine, parses the Tecplot file to fill a vtkMultiBlockDataSet object...
void Init()
This function initializes the context.
void GetArraysFromPointPackingZone(int numNodes, vtkPoints *theNodes, vtkPointData *nodeData)
This function extracts each variable array from a point-packing (tuple- based) zone and collects the ...
int GetNumberOfDataArrays()
Get the number of all data attributes (point data and cell data).
const char * GetDataTitle()
Get the Tecplot data title.
int IsDataAttributeCellBased(const char *attrName)
Get the type (0 for node-based and 1 for cell-based) of a specified data attribute (not 3D coordinate...
void GetPolyhedralGridCells(int numberCells, int numFaces, vtkUnstructuredGrid *unstruct) const
This function fills an allocated vtkUnstructuredGrid object with numberCells polyhedral cells to defi...
void GetStructuredGridFromBlockPackingZone(int iDimSize, int jDimSize, int kDimSize, int zoneIndx, const char *zoneName, vtkMultiBlockDataSet *multZone)
This function creates a vtkStructuredGrid object made up of a set of points and the associated attrib...
int IsDataAttributeCellBased(int attrIndx)
Get the type (0 for node-based and 1 for cell-based) of a specified data attribute (not 3D coordinate...
void GetPolygonalGridFromBlockPackingZone(int numNodes, int numElements, int numFaces, int zoneIndex, const char *zoneName, vtkMultiBlockDataSet *multZone)
This function creates a polygonal vtkUnstructuredGrid object made up of a set of points and the assoc...
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkDataArraySelection * DataArraySelection
const char * GetDataAttributeName(int attrIndx)
Get the name of a zero-based data attribute (not 3D coordinates).
void GetPolyhedralGridFromBlockPackingZone(int numNodes, int numElements, int numFaces, int zoneIndex, const char *zoneName, vtkMultiBlockDataSet *multZone)
This function creates a polyhedral vtkUnstructuredGrid object made up of a set of points and the asso...
void SetFileName(VTK_FILEPATH const char *fileName)
Specify a Tecplot ASCII file for data loading.
int GetNumberOfDataAttributes()
Get the number of standard data attributes (node-based and cell-based), excluding 3D coordinates.
vtkTecplotReaderInternal * Internal
std::string DataTitle
void GetUnstructuredGridFromBlockPackingZone(int numNodes, int numCells, const char *cellType, int zoneIndx, const char *zoneName, vtkMultiBlockDataSet *multZone)
This function creates a vtkUnstructuredGrid object made up of a set of points and the associated attr...
const char * GetDataArrayName(int arrayIdx)
Get the name of a data array specified by the zero-based index (arrayIdx).
vtkCallbackCommand * SelectionObserver
void GetPolygonalGridCells(int numFaces, int numEdges, vtkUnstructuredGrid *unstruct) const
This function fills an allocated vtkUnstructuredGrid object with numberCells polygonal cells to defin...
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
~vtkTecplotReader() override
std::vector< std::string > Variables
std::vector< std::string > ZoneNames
static void SelectionModifiedCallback(vtkObject *, unsigned long, void *tpReader, void *)
A callback function registered with the selection observer.
void SetDataArrayStatus(const char *arayName, int bChecked)
Set the status of a specific data array (0: de-select; 1: select) specified by the name.
std::vector< int > CellBased
void GetArraysFromBlockPackingZone(int numNodes, int numCells, vtkPoints *theNodes, vtkPointData *nodeData, vtkCellData *cellData)
This function extracts each variable array from a block-packing (component- based) zone and collects ...
const char * GetBlockName(int blockIdx)
Get the name of a block specified by a zero-based index.
static vtkTecplotReader * New()
void GetUnstructuredGridFromPointPackingZone(int numNodes, int numCells, const char *cellType, int zoneIndx, const char *zoneName, vtkMultiBlockDataSet *multZone)
This function creates a vtkUnstructuredGrid object made up of a set of points and the associated attr...
void GetStructuredGridFromPointPackingZone(int iDimSize, int jDimSize, int kDimSize, int zoneIndx, const char *zoneName, vtkMultiBlockDataSet *multZone)
This function creates a vtkStructuredGrid object made up of a set of points and the associated attrib...
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
int GetNumberOfBlocks()
Get the number of blocks (i.e., zones in Tecplot terms).
void GetUnstructuredGridCells(int numberCells, const char *cellTypeStr, vtkUnstructuredGrid *unstrctGrid)
This function fills an allocated vtkUnstructuredGrid object with numberCells cells of type cellTypeSt...
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
void GetDataArraysList()
Get the data arrays list from the tecplot file header.
int GetDataArrayStatus(const char *arayName)
Get the status of a specific data array (0: un-selected; 1: selected).
dataset represents arbitrary combinations of all possible cell types
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
@ string
Definition: vtkX3D.h:490
#define VTK_FILEPATH