VTK
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 
79 #ifndef vtkTecplotReader_h
80 #define vtkTecplotReader_h
81 
82 #include "vtkIOGeometryModule.h" // For export macro
84 
85 #include <vector> // STL Header; Required for vector
86 #include <string> // STL Header; Required for string
87 
88 class vtkPoints;
89 class vtkCellData;
90 class vtkPointData;
91 class vtkCallbackCommand;
95 class vtkTecplotReaderInternal;
96 
97 class VTKIOGEOMETRY_EXPORT vtkTecplotReader : public vtkMultiBlockDataSetAlgorithm
98 {
99 public:
100  static vtkTecplotReader * New();
102  void PrintSelf( ostream & os, vtkIndent indent );
103 
105 
108  vtkGetMacro( NumberOfVariables, int );
110 
114  void SetFileName( const char * fileName );
115 
119  const char * GetDataTitle();
120 
125 
130  const char * GetBlockName( int blockIdx );
131 
137 
142  const char * GetDataAttributeName( int attrIndx );
143 
149  int IsDataAttributeCellBased( const char * attrName );
150 
156  int IsDataAttributeCellBased( int attrIndx );
157 
162 
166  const char * GetDataArrayName( int arrayIdx );
167 
171  int GetDataArrayStatus( const char * arayName );
172 
177  void SetDataArrayStatus( const char * arayName, int bChecked );
178 
179 protected:
182 
184  virtual int RequestInformation( vtkInformation * request,
185  vtkInformationVector ** inputVector,
186  vtkInformationVector * outputVector );
187  virtual int RequestData
189 
194  ( vtkObject *, unsigned long, void * tpReader, void * );
195 
201  void Init();
202 
207 
212  void ReadFile( vtkMultiBlockDataSet * multZone);
213 
220  void GetArraysFromBlockPackingZone( int numNodes, int numCells,
221  vtkPoints * theNodes, vtkPointData * nodeData, vtkCellData * cellData );
222 
232  ( int numNodes, vtkPoints * theNodes, vtkPointData * nodeData );
233 
241  void GetStructuredGridFromBlockPackingZone( int iDimSize, int jDimSize,
242  int kDimSize, int zoneIndx, const char * zoneName,
243  vtkMultiBlockDataSet * multZone );
244 
252  void GetStructuredGridFromPointPackingZone( int iDimSize, int jDimSize,
253  int kDimSize, int zoneIndx, const char * zoneName,
254  vtkMultiBlockDataSet * multZone );
255 
263  void GetUnstructuredGridFromBlockPackingZone( int numNodes, int numCells,
264  const char * cellType, int zoneIndx, const char * zoneName,
265  vtkMultiBlockDataSet * multZone );
266 
274  void GetUnstructuredGridFromPointPackingZone( int numNodes, int numCells,
275  const char * cellType,int zoneIndx, const char * zoneName,
276  vtkMultiBlockDataSet * multZone );
277 
282  void GetUnstructuredGridCells( int numberCells, const char * cellTypeStr,
283  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 
298  vtkTecplotReader( const vtkTecplotReader & ) VTK_DELETE_FUNCTION;
299  void operator = ( const vtkTecplotReader & ) VTK_DELETE_FUNCTION;
300 };
301 
302 #endif
supports function callbacks
represent and manipulate cell attribute data
Definition: vtkCellData.h:39
Store on/off settings for data arrays for a vtkSource.
a simple class to control print indentation
Definition: vtkIndent.h:40
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:60
represent and manipulate point attribute data
Definition: vtkPointData.h:38
represent and manipulate 3D points
Definition: vtkPoints.h:40
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...
virtual int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
void Init()
This function initializes the context.
void SetFileName(const char *fileName)
Specify a Tecplot ASCII file for data loading.
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).
std::vector< std::string > ZoneNames
void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
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 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...
vtkDataArraySelection * DataArraySelection
const char * GetDataAttributeName(int attrIndx)
Get the name of a zero-based data attribute (not 3D coordinates).
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
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
std::vector< std::string > Variables
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.
virtual int FillOutputPortInformation(int port, vtkInformation *info)
Fill the output port information objects for this algorithm.
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...
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
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 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...
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
CellTypeInDataSet cellType(vtkDataSet *input)
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
@ string
Definition: vtkX3D.h:490