VTK  9.1.0
vtkTecplotReader.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkTecplotReader.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=========================================================================*/
15
16/*****************************************************************************
17 *
18 * Copyright (c) 2000 - 2009, Lawrence Livermore National Security, LLC
19 * Produced at the Lawrence Livermore National Laboratory
20 * LLNL-CODE-400124
21 * All rights reserved.
22 *
23 * This file was adapted from the ASCII Tecplot reader of VisIt. For details,
24 * see https://visit.llnl.gov/. The full copyright notice is contained in the
25 * file COPYRIGHT located at the root of the VisIt distribution or at
26 * http://www.llnl.gov/visit/copyright.html.
27 *
28 *****************************************************************************/
29
76#ifndef vtkTecplotReader_h
77#define vtkTecplotReader_h
78
79#include "vtkIOGeometryModule.h" // For export macro
81
82#include <string> // STL Header; Required for string
83#include <vector> // STL Header; Required for vector
84
85class vtkPoints;
86class vtkCellData;
87class vtkPointData;
92class vtkTecplotReaderInternal;
93
94class VTKIOGEOMETRY_EXPORT vtkTecplotReader : public vtkMultiBlockDataSetAlgorithm
95{
96public:
99 void PrintSelf(ostream& os, vtkIndent indent) override;
100
102
105 vtkGetMacro(NumberOfVariables, int);
107
111 void SetFileName(VTK_FILEPATH const char* fileName);
112
116 const char* GetDataTitle();
117
122
127 const char* GetBlockName(int blockIdx);
128
134
139 const char* GetDataAttributeName(int attrIndx);
140
146 int IsDataAttributeCellBased(const char* attrName);
147
153 int IsDataAttributeCellBased(int attrIndx);
154
159
163 const char* GetDataArrayName(int arrayIdx);
164
168 int GetDataArrayStatus(const char* arayName);
169
174 void SetDataArrayStatus(const char* arayName, int bChecked);
175
176protected:
179
182 vtkInformationVector* outputVector) override;
184
188 static void SelectionModifiedCallback(vtkObject*, unsigned long, void* tpReader, void*);
189
195 void Init();
196
201
207
215 int numNodes, int numCells, vtkPoints* theNodes, vtkPointData* nodeData, vtkCellData* cellData);
216
225 void GetArraysFromPointPackingZone(int numNodes, vtkPoints* theNodes, vtkPointData* nodeData);
226
234 void GetStructuredGridFromBlockPackingZone(int iDimSize, int jDimSize, int kDimSize, int zoneIndx,
235 const char* zoneName, vtkMultiBlockDataSet* multZone);
236
244 void GetStructuredGridFromPointPackingZone(int iDimSize, int jDimSize, int kDimSize, int zoneIndx,
245 const char* zoneName, vtkMultiBlockDataSet* multZone);
246
254 void GetUnstructuredGridFromBlockPackingZone(int numNodes, int numCells, const char* cellType,
255 int zoneIndx, const char* zoneName, vtkMultiBlockDataSet* multZone);
256
264 void GetPolyhedralGridFromBlockPackingZone(int numNodes, int numElements, int numFaces,
265 int zoneIndex, const char* zoneName, vtkMultiBlockDataSet* multZone);
266
274 void GetPolygonalGridFromBlockPackingZone(int numNodes, int numElements, int numFaces,
275 int zoneIndex, const char* zoneName, vtkMultiBlockDataSet* multZone);
276
281 void GetPolyhedralGridCells(int numberCells, int numFaces, vtkUnstructuredGrid* unstruct) const;
282
287 void GetPolygonalGridCells(int numFaces, int numEdges, vtkUnstructuredGrid* unstruct) const;
288
296 void GetUnstructuredGridFromPointPackingZone(int numNodes, int numCells, const char* cellType,
297 int zoneIndx, const char* zoneName, vtkMultiBlockDataSet* multZone);
298
304 int numberCells, const char* cellTypeStr, vtkUnstructuredGrid* unstrctGrid);
305
307 char* FileName;
310 vtkTecplotReaderInternal* Internal;
311
313 std::vector<int> CellBased;
314 std::vector<std::string> ZoneNames;
315 std::vector<std::string> Variables;
316
317private:
318 vtkTecplotReader(const vtkTecplotReader&) = delete;
319 void operator=(const vtkTecplotReader&) = delete;
320};
321
322#endif
supports function callbacks
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
Store on/off settings for data arrays for a vtkSource.
a simple class to control print indentation
Definition: vtkIndent.h:113
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:73
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
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).
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.
const char * GetDataTitle()
Get the Tecplot data title.
vtkDataArraySelection * DataArraySelection
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...
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
static vtkTecplotReader * New()
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 ...
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).
const char * GetDataArrayName(int arrayIdx)
Get the name of a data array specified by the zero-based index (arrayIdx).
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.
const char * GetBlockName(int blockIdx)
Get the name of a block specified by a zero-based index.
void GetDataArraysList()
Get the data arrays list from the tecplot file header.
const char * GetDataAttributeName(int attrIndx)
Get the name of a zero-based data attribute (not 3D coordinates).
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:382
@ port
Definition: vtkX3D.h:453
@ string
Definition: vtkX3D.h:496
#define VTK_FILEPATH