VTK
|
00001 /*========================================================================= 00002 00003 Program: Visualization Toolkit 00004 Module: vtkExodusReader.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 Copyright (c) Sandia Corporation 00017 See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details. 00018 ----------------------------------------------------------------------------*/ 00019 00039 #ifndef __vtkExodusReader_h 00040 #define __vtkExodusReader_h 00041 00042 #define ARRAY_TYPE_NAMES_IN_CXX_FILE 00043 00044 #include "vtkUnstructuredGridAlgorithm.h" 00045 00046 class vtkIntArray; 00047 class vtkFloatArray; 00048 class vtkDataArray; 00049 class vtkDataSet; 00050 class vtkPoints; 00051 class vtkExodusMetadata; 00052 class vtkExodusModel; 00053 class vtkExodusXMLParser; 00054 00055 00056 #include "vtkDSPFilterGroup.h" //for USE_EXO_DSP_FILTERS 00057 00058 00059 class VTK_HYBRID_EXPORT vtkExodusReader : public vtkUnstructuredGridAlgorithm 00060 { 00061 public: 00062 static vtkExodusReader *New(); 00063 vtkTypeMacro(vtkExodusReader,vtkUnstructuredGridAlgorithm); 00064 void PrintSelf(ostream& os, vtkIndent indent); 00065 00067 int CanReadFile(const char* fname); 00068 00070 00072 VTK_LEGACY(void SetFileName(const char* fname)); 00073 vtkGetStringMacro(FileName); 00075 00077 00078 vtkSetStringMacro(XMLFileName); 00079 vtkGetStringMacro(XMLFileName); 00081 00083 00084 vtkSetMacro(TimeStep, int); 00085 vtkGetMacro(TimeStep, int); 00087 00089 00092 vtkSetMacro(GenerateBlockIdCellArray, int); 00093 vtkGetMacro(GenerateBlockIdCellArray, int); 00094 vtkBooleanMacro(GenerateBlockIdCellArray, int); 00095 const char *GetBlockIdArrayName() { return "BlockId"; } 00097 00098 00100 00103 vtkSetMacro(GenerateGlobalElementIdArray, int); 00104 vtkGetMacro(GenerateGlobalElementIdArray, int); 00105 vtkBooleanMacro(GenerateGlobalElementIdArray, int); 00106 //BTX 00107 enum { 00108 SEARCH_TYPE_ELEMENT=0, 00109 SEARCH_TYPE_NODE, 00110 SEARCH_TYPE_ELEMENT_THEN_NODE, 00111 SEARCH_TYPE_NODE_THEN_ELEMENT, 00112 ID_NOT_FOUND=-234121312 00113 }; 00114 //ETX 00115 static const char *GetGlobalElementIdArrayName() { return "GlobalElementId"; } 00116 static const char *GetPedigreeElementIdArrayName() { return "PedigreeElementId"; } 00117 static int GetGlobalElementID( vtkDataSet *data, int localID ); 00118 static int GetGlobalElementID ( vtkDataSet *data, int localID, 00119 int searchType ); 00121 00123 00127 vtkSetMacro(GenerateGlobalNodeIdArray, int); 00128 vtkGetMacro(GenerateGlobalNodeIdArray, int); 00129 vtkBooleanMacro(GenerateGlobalNodeIdArray, int); 00130 static const char *GetGlobalNodeIdArrayName() { return "GlobalNodeId"; } 00131 static const char *GetPedigreeNodeIdArrayName() { return "PedigreeNodeId"; } 00132 static int GetGlobalNodeID( vtkDataSet *data, int localID ); 00133 static int GetGlobalNodeID( vtkDataSet *data, int localID, 00134 int searchType ); 00136 00138 00142 vtkSetMacro(ApplyDisplacements, int); 00143 vtkGetMacro(ApplyDisplacements, int); 00144 vtkBooleanMacro(ApplyDisplacements, int); 00145 vtkSetMacro(DisplacementMagnitude, float); 00146 vtkGetMacro(DisplacementMagnitude, float); 00148 00150 00151 vtkGetStringMacro(Title); 00152 vtkGetMacro(Dimensionality, int); 00153 vtkGetMacro(NumberOfTimeSteps, int); 00154 int GetNumberOfElements() { return this->NumberOfUsedElements; } 00155 vtkGetMacro(NumberOfNodeSets, int); 00156 vtkGetMacro(NumberOfSideSets, int); 00157 vtkGetMacro(NumberOfBlocks, int); 00158 vtkGetVector2Macro(TimeStepRange, int); 00159 vtkSetVector2Macro(TimeStepRange, int); 00160 int GetNumberOfNodes() { return this->NumberOfUsedNodes; } 00161 int GetNumberOfElementsInBlock(int block_idx); 00162 int GetBlockId(int block_idx); 00163 virtual int GetTotalNumberOfNodes() { return this->NumberOfNodesInFile; } 00165 00166 00168 00173 int GetNumberOfPointArrays(); 00174 const char *GetPointArrayName(int index); 00175 int GetPointArrayID( const char *name ); 00176 int GetPointArrayNumberOfComponents(int index); 00177 void SetPointArrayStatus(int index, int flag); 00178 void SetPointArrayStatus(const char*, int flag); 00179 int GetPointArrayStatus(int index); 00180 int GetPointArrayStatus(const char*); 00182 00183 int GetNumberOfCellArrays(); 00184 const char *GetCellArrayName(int index); 00185 int GetCellArrayID( const char *name ); 00186 int GetCellArrayNumberOfComponents(int index); 00187 void SetCellArrayStatus(int index, int flag); 00188 void SetCellArrayStatus(const char*, int flag); 00189 int GetCellArrayStatus(int index); 00190 int GetCellArrayStatus(const char*); 00191 virtual int GetTotalNumberOfElements() 00192 { return this->NumberOfElementsInFile; } 00193 00195 00199 int GetNumberOfBlockArrays(); 00200 const char *GetBlockArrayName(int index); 00201 int GetBlockArrayID( const char *name ); 00202 void SetBlockArrayStatus(int index, int flag); 00203 void SetBlockArrayStatus(const char*, int flag); 00204 int GetBlockArrayStatus(int index); 00205 int GetBlockArrayStatus(const char*); 00207 00208 00210 00217 int GetNumberOfNodeSetArrays(){return this->GetNumberOfNodeSets();} 00218 int GetNodeSetArrayStatus(int index); 00219 int GetNodeSetArrayStatus(const char* name); 00220 void SetNodeSetArrayStatus(int index, int flag); 00221 void SetNodeSetArrayStatus(const char* name, int flag); 00222 const char *GetNodeSetArrayName(int index); 00224 00225 int GetNumberOfSideSetArrays(){return this->GetNumberOfSideSets();} 00226 int GetSideSetArrayStatus(int index); 00227 int GetSideSetArrayStatus(const char* name); 00228 void SetSideSetArrayStatus(int index, int flag); 00229 void SetSideSetArrayStatus(const char* name, int flag); 00230 const char *GetSideSetArrayName(int index); 00231 00233 00237 int GetNumberOfPartArrays(); 00238 const char *GetPartArrayName(int arrayIdx); 00239 int GetPartArrayID( const char *name ); 00240 const char *GetPartBlockInfo(int arrayIdx); 00241 void SetPartArrayStatus(int index, int flag); 00242 void SetPartArrayStatus(const char*, int flag); 00243 int GetPartArrayStatus(int index); 00244 int GetPartArrayStatus(const char*); 00246 00247 00249 00253 int GetNumberOfMaterialArrays(); 00254 const char *GetMaterialArrayName(int arrayIdx); 00255 int GetMaterialArrayID( const char *name ); 00256 void SetMaterialArrayStatus(int index, int flag); 00257 void SetMaterialArrayStatus(const char*, int flag); 00258 int GetMaterialArrayStatus(int index); 00259 int GetMaterialArrayStatus(const char*); 00261 00263 00267 int GetNumberOfAssemblyArrays(); 00268 const char *GetAssemblyArrayName(int arrayIdx); 00269 int GetAssemblyArrayID( const char *name ); 00270 void SetAssemblyArrayStatus(int index, int flag); 00271 void SetAssemblyArrayStatus(const char*, int flag); 00272 int GetAssemblyArrayStatus(int index); 00273 int GetAssemblyArrayStatus(const char*); 00275 00277 00284 int GetNumberOfHierarchyArrays(); 00285 const char *GetHierarchyArrayName(int arrayIdx); 00286 void SetHierarchyArrayStatus(int index, int flag); 00287 void SetHierarchyArrayStatus(const char*, int flag); 00288 int GetHierarchyArrayStatus(int index); 00289 int GetHierarchyArrayStatus(const char*); 00291 00293 00300 vtkGetMacro(HasModeShapes, int); 00301 vtkSetMacro(HasModeShapes, int); 00302 vtkBooleanMacro(HasModeShapes, int); 00304 00305 vtkGetMacro(DisplayType,int); 00306 virtual void SetDisplayType(int type); 00307 00313 vtkBooleanMacro(ExodusModelMetadata, int); 00314 vtkSetMacro(ExodusModelMetadata, int); 00315 vtkGetMacro(ExodusModelMetadata, int); 00316 00319 vtkExodusModel *GetExodusModel(){return this->ExodusModel;} 00320 00328 vtkSetMacro(PackExodusModelOntoOutput, int); 00329 vtkGetMacro(PackExodusModelOntoOutput, int); 00330 vtkBooleanMacro(PackExodusModelOntoOutput, int); 00331 00332 //BTX 00334 00335 enum ArrayType { 00336 CELL=0, 00337 POINT, 00338 BLOCK, 00339 PART, 00340 MATERIAL, 00341 ASSEMBLY, 00342 HIERARCHY, 00343 NUM_ARRAY_TYPES, 00344 UNKNOWN_TYPE 00345 }; 00346 //ETX 00348 00350 int IsValidVariable( const char *type, const char *name ); 00351 00352 //BTX 00354 00355 int GetNumberOfArrays( vtkExodusReader::ArrayType type ); 00356 const char *GetArrayName( vtkExodusReader::ArrayType type, int id ); 00357 //ETX 00359 00361 int GetVariableID ( const char *type, const char *name ); 00362 00363 void SetAllAssemblyArrayStatus( int status ); 00364 void SetAllBlockArrayStatus( int status ); 00365 void SetAllCellArrayStatus( int status ); 00366 void SetAllHierarchyArrayStatus( int status ); 00367 void SetAllMaterialArrayStatus( int status ); 00368 void SetAllPartArrayStatus( int status ); 00369 void SetAllPointArrayStatus( int status ); 00370 //BTX 00371 void SetAllArrayStatus ( vtkExodusReader::ArrayType type, int flag ); 00372 void SetArrayStatus ( vtkExodusReader::ArrayType type, const char *name, 00373 int flag ); 00374 //ETX 00375 void SetArrayStatus ( const char *type, const char *name, int flag ) 00376 { 00377 this->SetArrayStatus( this->GetArrayTypeID(type), name, flag ); 00378 } 00379 //BTX 00380 int GetArrayStatus ( vtkExodusReader::ArrayType type, const char *name ); 00381 //ETX 00382 int GetArrayStatus ( const char *type, const char *name ) 00383 { 00384 return this->GetArrayStatus( this->GetArrayTypeID( type ), name ); 00385 } 00386 00387 // Helper functions 00388 static int StringsEqual(const char* s1, char* s2); 00389 static void StringUppercase(const char* str, char* upperstr); 00390 static char *StrDupWithNew(const char *s); 00391 00392 // time series query functions 00393 int GetTimeSeriesData( int ID, const char *vName, const char *vType, 00394 vtkFloatArray *result ); 00395 00396 00397 //begin USE_EXO_DSP_FILTERS 00398 int GetNumberOfVariableArrays(); 00399 const char *GetVariableArrayName(int a_which); 00400 void EnableDSPFiltering(); 00401 void AddFilter(vtkDSPFilterDefinition *a_filter); 00402 void StartAddingFilter(); 00403 void AddFilterInputVar(char *name); 00404 void AddFilterOutputVar(char *name); 00405 void AddFilterNumeratorWeight(double weight); 00406 void AddFilterForwardNumeratorWeight(double weight); 00407 void AddFilterDenominatorWeight(double weight); 00408 void FinishAddingFilter(); 00409 void RemoveFilter(char *a_outputVariableName); 00410 void GetDSPOutputArrays(int exoid, vtkUnstructuredGrid* output); 00411 //BTX 00412 vtkExodusReader::ArrayType GetArrayTypeID( const char *type ); 00413 00414 #ifdef ARRAY_TYPE_NAMES_IN_CXX_FILE 00415 static const char *GetArrayTypeName( vtkExodusReader::ArrayType type ); 00416 #else 00417 static const char *ArrayTypeNames[NUM_ARRAY_TYPES]; 00418 00419 static const char *GetArrayTypeName( vtkExodusReader::ArrayType type ) 00420 { 00421 return ArrayTypeNames[type]; 00422 } 00423 #endif 00424 //ETX 00425 00426 vtkDSPFilterDefinition *AddingFilter; 00427 int DSPFilteringIsEnabled; 00428 vtkDSPFilterGroup **DSPFilters; 00429 //end USE_EXO_DSP_FILTERS 00430 00431 00432 00433 protected: 00434 vtkExodusReader(); 00435 ~vtkExodusReader(); 00436 00437 void NewExodusModel(); 00438 00439 void ReadGeometry(int exoid, vtkUnstructuredGrid* output); 00440 void ReadCells(int exoid, vtkUnstructuredGrid* output); 00441 void ReadPoints(int exoid, vtkUnstructuredGrid* output); 00442 void ReadArrays(int exoid, vtkUnstructuredGrid* output); 00443 void ReadNodeAndSideSets(int exoid, vtkUnstructuredGrid* output); 00444 vtkDataArray *ReadPointArray(int exoid, int varIndex); 00445 vtkDataArray *ReadPointVector(int handle, int varIndex, int dim); 00446 vtkDataArray *ReadCellArray(int exoid, int varIndex); 00447 vtkDataArray *ReadCellVector(int handle, int varIndex, int dim); 00448 void ReadNodeSetMetadata(); 00449 void ReadSideSetMetadata(); 00450 00451 // helper for finding IDs 00452 static int GetIDHelper ( const char *arrayName, vtkDataSet *data, int localID, 00453 int searchType ); 00454 static int GetGlobalID( const char *arrayName, vtkDataSet *data, int localID, 00455 int searchType ); 00456 00457 // This method is a helper for determining the 00458 // number of additional cell scalar field 00459 // values needed to 'pad' for node and side sets 00460 int GetExtraCellCountForNodeSideSets(); 00461 00462 // This method generates arrays like blockid, global nodeid 00463 // and global element id 00464 void GenerateExtraArrays(vtkUnstructuredGrid* output); 00465 00466 // Parameters for controlling what is read in. 00467 char *FileName; 00468 char *XMLFileName; 00469 int TimeStep; 00470 int ActualTimeStep; 00471 double TimeValue; 00472 int GenerateBlockIdCellArray; 00473 int GenerateGlobalElementIdArray; 00474 int GenerateGlobalNodeIdArray; 00475 int ApplyDisplacements; 00476 double DisplacementMagnitude; 00477 00478 // Information specific for exodus files. 00479 vtkSetStringMacro(Title); 00480 char *Title; 00481 int Dimensionality; 00482 int NumberOfNodeSets; 00483 int NumberOfSideSets; 00484 int NumberOfBlocks; 00485 int NumberOfUsedNodes; 00486 int NumberOfNodesInFile; 00487 int NumberOfUsedElements; 00488 int NumberOfElementsInFile; 00489 int NumberOfTimeSteps; 00490 int ExodusCPUWordSize; 00491 int ExodusIOWordSize; 00492 float ExodusVersion; 00493 vtkIntArray *CellVarTruthTable; 00494 00495 //1=display Block names 00496 //2=display Part names 00497 //3=display Material names 00498 int DisplayType; 00499 00500 //Parser that understands the xml part and material file 00501 vtkExodusXMLParser *Parser; 00502 00503 // **KEN** By VTK convention, metaData should be Metadata. 00504 00506 // Scalar Array and Block Info 00508 vtkExodusMetadata *MetaData; 00509 00510 00512 00513 int CurrentHandle; 00514 char* CurrentFileName; 00515 char* CurrentXMLFileName; 00516 vtkSetStringMacro(CurrentFileName); 00517 vtkSetStringMacro(CurrentXMLFileName); 00519 00520 // Open the exodus file, and set some basic information 00521 int OpenCurrentFile(); 00522 00523 // Close the exodus file 00524 void CloseCurrentFile(); 00525 00526 00528 int TimeStepRange[2]; 00529 00530 // DataCache: this object keeps the points and cells 00531 // around so they don't need to be re-read when the 00532 // timestep changes or an scalar array is switched 00533 vtkUnstructuredGrid *DataCache; 00534 00535 // Should I re-read in the geometry and topology of the dataset 00536 int RemakeDataCacheFlag; 00537 00538 // vtkExodusModel needs to count changes in geometry, so it knows 00539 // if geometry has changed since it last updated model data. 00540 00541 int NewGeometryCount; 00542 00543 // PointMap keeps track of which points are actually 00544 // used by the cells that are read in (blocks) 00545 vtkIntArray *PointMap; 00546 vtkIntArray *ReversePointMap; 00547 void SetUpPointMap(int num_points); 00548 int GetPointMapIndex(int point_id); 00549 00550 // Global element ID cache 00551 int *GlobalElementIdCache; 00552 void SetGlobalElementIdCache(int *list); 00553 00554 // Time query function. Called by ExecuteInformation(). 00555 // Fills the TimestepValues array. 00556 void GetAllTimes(vtkInformationVector *); 00557 00558 int HasModeShapes; 00559 00560 vtkExodusModel *ExodusModel; 00561 int PackExodusModelOntoOutput; 00562 int ExodusModelMetadata; 00563 00564 double *TimeSteps; 00565 00566 int RequestInformation( 00567 vtkInformation *, vtkInformationVector **, vtkInformationVector *); 00568 int RequestData( 00569 vtkInformation *, vtkInformationVector **, vtkInformationVector *); 00570 00571 // Used to determine current progress. 00572 double ProgressOffset; 00573 double ProgressScale; 00574 00575 private: 00576 vtkExodusReader(const vtkExodusReader&); // Not implemented 00577 void operator=(const vtkExodusReader&); // Not implemented 00578 00579 void AddDisplacements(vtkUnstructuredGrid* output); 00580 void RemoveBeginningAndTrailingSpaces(char **names, int len); 00581 00582 void FixMetadataTruthTable(int *table, int len); 00583 }; 00584 00585 #endif