VTK  9.0.1
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 
85 class vtkPoints;
86 class vtkCellData;
87 class vtkPointData;
88 class vtkCallbackCommand;
92 class vtkTecplotReaderInternal;
93 
94 class VTKIOGEOMETRY_EXPORT vtkTecplotReader : public vtkMultiBlockDataSetAlgorithm
95 {
96 public:
97  static vtkTecplotReader* New();
99  void PrintSelf(ostream& os, vtkIndent indent) override;
100 
102 
105  vtkGetMacro(NumberOfVariables, int);
107 
111  void SetFileName(const char* fileName);
112 
116  const char* GetDataTitle();
117 
121  int GetNumberOfBlocks();
122 
127  const char* GetBlockName(int blockIdx);
128 
133  int GetNumberOfDataAttributes();
134 
139  const char* GetDataAttributeName(int attrIndx);
140 
146  int IsDataAttributeCellBased(const char* attrName);
147 
153  int IsDataAttributeCellBased(int attrIndx);
154 
158  int GetNumberOfDataArrays();
159 
163  const char* GetDataArrayName(int arrayIdx);
164 
168  int GetDataArrayStatus(const char* arayName);
169 
174  void SetDataArrayStatus(const char* arayName, int bChecked);
175 
176 protected:
178  ~vtkTecplotReader() override;
179 
181  int RequestInformation(vtkInformation* request, vtkInformationVector** inputVector,
182  vtkInformationVector* outputVector) override;
184 
188  static void SelectionModifiedCallback(vtkObject*, unsigned long, void* tpReader, void*);
189 
195  void Init();
196 
200  void GetDataArraysList();
201 
206  void ReadFile(vtkMultiBlockDataSet* multZone);
207 
214  void GetArraysFromBlockPackingZone(
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 
303  void GetUnstructuredGridCells(
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 
317 private:
318  vtkTecplotReader(const vtkTecplotReader&) = delete;
319  void operator=(const vtkTecplotReader&) = delete;
320 };
321 
322 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:33
vtkTecplotReader
A concrete class to read an ASCII Tecplot file.
Definition: vtkTecplotReader.h:94
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:31
vtkInformationVector
Store zero or more vtkInformation instances.
Definition: vtkInformationVector.h:35
vtkTecplotReader::FileName
char * FileName
Definition: vtkTecplotReader.h:307
vtkMultiBlockDataSetAlgorithm::FillOutputPortInformation
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
vtkObject
abstract base class for most VTK objects
Definition: vtkObject.h:62
vtkMultiBlockDataSetAlgorithm::New
static vtkMultiBlockDataSetAlgorithm * New()
vtkTecplotReader::ZoneNames
std::vector< std::string > ZoneNames
Definition: vtkTecplotReader.h:314
vtkMultiBlockDataSet
Composite dataset that organizes datasets into blocks.
Definition: vtkMultiBlockDataSet.h:45
vtkTecplotReader::SelectionObserver
vtkCallbackCommand * SelectionObserver
Definition: vtkTecplotReader.h:308
vtkMultiBlockDataSetAlgorithm::RequestData
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
Definition: vtkMultiBlockDataSetAlgorithm.h:89
vtkDataArraySelection
Store on/off settings for data arrays for a vtkSource.
Definition: vtkDataArraySelection.h:34
vtkTecplotReader::Variables
std::vector< std::string > Variables
Definition: vtkTecplotReader.h:315
vtkX3D::port
@ port
Definition: vtkX3D.h:453
vtkTecplotReader::DataArraySelection
vtkDataArraySelection * DataArraySelection
Definition: vtkTecplotReader.h:309
vtkTecplotReader::NumberOfVariables
int NumberOfVariables
Definition: vtkTecplotReader.h:306
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkMultiBlockDataSetAlgorithm.h
vtkTecplotReader::CellBased
std::vector< int > CellBased
Definition: vtkTecplotReader.h:313
vtkMultiBlockDataSetAlgorithm::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkInformation
Store vtkAlgorithm input/output information.
Definition: vtkInformation.h:73
vtkX3D::info
@ info
Definition: vtkX3D.h:382
vtkX3D::string
@ string
Definition: vtkX3D.h:496
vtkCallbackCommand
supports function callbacks
Definition: vtkCallbackCommand.h:44
vtkTecplotReader::Internal
vtkTecplotReaderInternal * Internal
Definition: vtkTecplotReader.h:310
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:92
vtkMultiBlockDataSetAlgorithm::RequestInformation
virtual int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
Definition: vtkMultiBlockDataSetAlgorithm.h:80
vtkTecplotReader::DataTitle
std::string DataTitle
Definition: vtkTecplotReader.h:312
vtkMultiBlockDataSetAlgorithm
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
Definition: vtkMultiBlockDataSetAlgorithm.h:32