VTK
|
00001 /*========================================================================= 00002 00003 Program: Visualization Toolkit 00004 Module: vtkTecplotReader.h 00005 00006 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen 00007 All rights reserved. 00008 See Copyright.txt or http://www.kitware.com/Copyright.htm for details. 00009 00010 This software is distributed WITHOUT ANY WARRANTY; without even 00011 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 00012 PURPOSE. See the above copyright notice for more information. 00013 00014 =========================================================================*/ 00015 00016 /***************************************************************************** 00017 * 00018 * Copyright (c) 2000 - 2009, Lawrence Livermore National Security, LLC 00019 * Produced at the Lawrence Livermore National Laboratory 00020 * LLNL-CODE-400124 00021 * All rights reserved. 00022 * 00023 * This file was adapted from the ASCII Tecplot reader of VisIt. For details, 00024 * see https://visit.llnl.gov/. The full copyright notice is contained in the 00025 * file COPYRIGHT located at the root of the VisIt distribution or at 00026 * http://www.llnl.gov/visit/copyright.html. 00027 * 00028 *****************************************************************************/ 00029 00078 #ifndef __vtkTecplotReader_h 00079 #define __vtkTecplotReader_h 00080 00081 #include "vtkIOGeometryModule.h" // For export macro 00082 #include "vtkMultiBlockDataSetAlgorithm.h" 00083 00084 //BTX 00085 #include <vector> // STL Header; Required for vector 00086 #include <string> // STL Header; Required for string 00087 //ETX 00088 00089 class vtkPoints; 00090 class vtkCellData; 00091 class vtkPointData; 00092 class vtkCallbackCommand; 00093 class vtkUnstructuredGrid; 00094 class vtkMultiBlockDataSet; 00095 class vtkDataArraySelection; 00096 class vtkTecplotReaderInternal; 00097 00098 class VTKIOGEOMETRY_EXPORT vtkTecplotReader : public vtkMultiBlockDataSetAlgorithm 00099 { 00100 public: 00101 static vtkTecplotReader * New(); 00102 vtkTypeMacro( vtkTecplotReader, vtkMultiBlockDataSetAlgorithm ); 00103 void PrintSelf( ostream & os, vtkIndent indent ); 00104 00106 00107 vtkGetMacro( NumberOfVariables, int ); 00109 00111 void SetFileName( const char * fileName ); 00112 00114 const char * GetDataTitle(); 00115 00117 int GetNumberOfBlocks(); 00118 00121 const char * GetBlockName( int blockIdx ); 00122 00125 int GetNumberOfDataAttributes(); 00126 00129 const char * GetDataAttributeName( int attrIndx ); 00130 00134 int IsDataAttributeCellBased( const char * attrName ); 00135 00139 int IsDataAttributeCellBased( int attrIndx ); 00140 00142 int GetNumberOfDataArrays(); 00143 00146 const char * GetDataArrayName( int arrayIdx ); 00147 00150 int GetDataArrayStatus( const char * arayName ); 00151 00154 void SetDataArrayStatus( const char * arayName, int bChecked ); 00155 00156 protected: 00157 vtkTecplotReader(); 00158 ~vtkTecplotReader(); 00159 00160 virtual int FillOutputPortInformation( int port, vtkInformation * info ); 00161 virtual int RequestInformation( vtkInformation * request, 00162 vtkInformationVector ** inputVector, 00163 vtkInformationVector * outputVector ); 00164 virtual int RequestData 00165 ( vtkInformation *, vtkInformationVector **, vtkInformationVector * ); 00166 00168 00169 static void SelectionModifiedCallback 00170 ( vtkObject *, unsigned long, void * tpReader, void * ); 00172 00177 void Init(); 00178 00180 void GetDataArraysList(); 00181 00184 void ReadFile( vtkMultiBlockDataSet * multZone); 00185 00187 00192 void GetArraysFromBlockPackingZone( int numNodes, int numCells, 00193 vtkPoints * theNodes, vtkPointData * nodeData, vtkCellData * cellData ); 00195 00197 00203 void GetArraysFromPointPackingZone 00204 ( int numNodes, vtkPoints * theNodes, vtkPointData * nodeData ); 00206 00208 00213 void GetStructuredGridFromBlockPackingZone( int iDimSize, int jDimSize, 00214 int kDimSize, int zoneIndx, const char * zoneName, 00215 vtkMultiBlockDataSet * multZone ); 00217 00219 00224 void GetStructuredGridFromPointPackingZone( int iDimSize, int jDimSize, 00225 int kDimSize, int zoneIndx, const char * zoneName, 00226 vtkMultiBlockDataSet * multZone ); 00228 00230 00235 void GetUnstructuredGridFromBlockPackingZone( int numNodes, int numCells, 00236 const char * cellType, int zoneIndx, const char * zoneName, 00237 vtkMultiBlockDataSet * multZone ); 00239 00241 00246 void GetUnstructuredGridFromPointPackingZone( int numNodes, int numCells, 00247 const char * cellType,int zoneIndx, const char * zoneName, 00248 vtkMultiBlockDataSet * multZone ); 00250 00252 00254 void GetUnstructuredGridCells( int numberCells, const char * cellTypeStr, 00255 vtkUnstructuredGrid * unstrctGrid ); 00257 00258 int NumberOfVariables; 00259 char * FileName; 00260 vtkCallbackCommand * SelectionObserver; 00261 vtkDataArraySelection * DataArraySelection; 00262 vtkTecplotReaderInternal * Internal; 00263 00264 //BTX 00265 std::string DataTitle; 00266 std::vector< int > CellBased; 00267 std::vector< std::string > ZoneNames; 00268 std::vector< std::string > Variables; 00269 //ETX 00270 00271 private: 00272 00273 vtkTecplotReader( const vtkTecplotReader & ); // Not implemented. 00274 void operator = ( const vtkTecplotReader & ); // Not implemented. 00275 }; 00276 00277 #endif