VTK
dox/IO/Exodus/vtkExodusIIReaderPrivate.h
Go to the documentation of this file.
00001 #ifndef __vtkExodusIIReaderPrivate_h
00002 #define __vtkExodusIIReaderPrivate_h
00003 
00004 // Do not include this file directly. It is only for use
00005 // from inside the ExodusII reader and its descendants.
00006 
00007 #include "vtkToolkits.h" // make sure VTK_USE_PARALLEL is properly set
00008 #include "vtkExodusIICache.h"
00009 #include "vtksys/RegularExpression.hxx"
00010 
00011 #include <map>
00012 #include <vector>
00013 
00014 #include "vtk_exodusII.h"
00015 #include "vtkIOExodusModule.h" // For export macro
00016 class vtkExodusIIReaderParser;
00017 class vtkMutableDirectedGraph;
00018 
00022 class VTKIOEXODUS_EXPORT vtkExodusIIReaderPrivate : public vtkObject
00023 {
00024 public:
00025   static vtkExodusIIReaderPrivate* New();
00026   void PrintData( ostream& os, vtkIndent indent );
00027   vtkTypeMacro(vtkExodusIIReaderPrivate,vtkObject);
00028   //virtual void Modified();
00029 
00031   int OpenFile( const char* filename );
00032 
00034   int CloseFile();
00035 
00037   int RequestInformation();
00038 
00040   vtkMutableDirectedGraph* GetSIL()
00041     { return this->SIL; }
00042 
00044   int RequestData( vtkIdType timeStep, vtkMultiBlockDataSet* output );
00045 
00050   int SetUpEmptyGrid( vtkMultiBlockDataSet* output );
00051 
00063   void Reset();
00064 
00069   void ResetSettings();
00070 
00072   void ResetCache();
00073 
00075   void SetCacheSize(double size);
00076 
00078   vtkGetMacro(CacheSize, double);
00079 
00084   int GetNumberOfTimeSteps() { return (int) this->Times.size(); }
00085 
00087   vtkGetMacro(TimeStep,int);
00088 
00090   vtkSetMacro(TimeStep,int);
00091 
00094   vtkGetMacro(SqueezePoints,int);
00095 
00098   void SetSqueezePoints( int sp );
00099 
00102   vtkBooleanMacro(SqueezePoints,int);
00103 
00105   int GetNumberOfNodes();
00106 
00111   int GetNumberOfObjectsOfType( int otype );
00112 
00123   int GetNumberOfObjectArraysOfType( int otype );
00124 
00129   const char* GetObjectName( int otype, int i );
00130 
00135   int GetObjectId( int otype, int i );
00136 
00143   int GetObjectSize( int otype, int i );
00144 
00149   int GetObjectStatus( int otype, int i );
00150 
00156   int GetUnsortedObjectStatus( int otype, int i );
00157 
00162   void SetObjectStatus( int otype, int i, int stat );
00163 
00169   void SetUnsortedObjectStatus( int otype, int i, int stat );
00170 
00175   const char* GetObjectArrayName( int otype, int i );
00176 
00181   int GetNumberOfObjectArrayComponents( int otype, int i );
00182 
00187   int GetObjectArrayStatus( int otype, int i );
00188 
00193   void SetObjectArrayStatus( int otype, int i, int stat );
00194 
00201   int GetNumberOfObjectAttributes( int objectType, int objectIndex );
00202   const char* GetObjectAttributeName( int objectType,
00203                                       int objectIndex,
00204                                       int attributeIndex );
00205   int GetObjectAttributeIndex( int objectType,
00206                                int objectIndex,
00207                                const char* attribName );
00208   int GetObjectAttributeStatus( int objectType,
00209                                 int objectIndex,
00210                                 int attribIndex );
00211   void SetObjectAttributeStatus( int objectType,
00212                                  int objectIndex,
00213                                  int attribIndex, int status );
00214 
00216   vtkGetMacro(GenerateObjectIdArray,int);
00217   vtkSetMacro(GenerateObjectIdArray,int);
00218   static const char* GetObjectIdArrayName() { return "ObjectId"; }
00219 
00220   vtkSetMacro(GenerateGlobalElementIdArray,int);
00221   vtkGetMacro(GenerateGlobalElementIdArray,int);
00222   static const char* GetGlobalElementIdArrayName() { return "GlobalElementId"; }
00223 
00224   vtkSetMacro(GenerateGlobalNodeIdArray,int);
00225   vtkGetMacro(GenerateGlobalNodeIdArray,int);
00226   static const char* GetGlobalNodeIdArrayName() { return "GlobalNodeId"; }
00227 
00228   vtkSetMacro(GenerateImplicitElementIdArray,int);
00229   vtkGetMacro(GenerateImplicitElementIdArray,int);
00230   static const char* GetImplicitElementIdArrayName() { return "ImplicitElementId"; }
00231 
00232   vtkSetMacro(GenerateImplicitNodeIdArray,int);
00233   vtkGetMacro(GenerateImplicitNodeIdArray,int);
00234   static const char* GetImplicitNodeIdArrayName() { return "ImplicitNodeId"; }
00235 
00239   vtkSetMacro(GenerateFileIdArray,int);
00240   vtkGetMacro(GenerateFileIdArray,int);
00241   static const char* GetFileIdArrayName() { return "FileId"; }
00242 
00244   vtkSetMacro(FileId,int);
00245   vtkGetMacro(FileId,int);
00246 
00247   static const char *GetGlobalVariableValuesArrayName()
00248     { return "GlobalVariableValues"; }
00249   static const char *GetGlobalVariableNamesArrayName()
00250     { return "GlobalVariableNames"; }
00251 
00252   virtual void SetApplyDisplacements( int d );
00253   vtkGetMacro(ApplyDisplacements,int);
00254 
00255   virtual void SetDisplacementMagnitude( double s );
00256   vtkGetMacro(DisplacementMagnitude,double);
00257 
00258   vtkSetMacro(HasModeShapes,int);
00259   vtkGetMacro(HasModeShapes,int);
00260 
00261   vtkSetMacro(ModeShapeTime,double);
00262   vtkGetMacro(ModeShapeTime,double);
00263 
00264   vtkSetMacro(AnimateModeShapes, int);
00265   vtkGetMacro(AnimateModeShapes, int);
00266 
00267   vtkDataArray* FindDisplacementVectors( int timeStep );
00268 
00269   const struct ex_init_params* GetModelParams() const
00270     { return &this->ModelParameters; }
00271 
00273   struct VTKIOEXODUS_EXPORT ArrayInfoType {
00275     vtkStdString Name;
00277     int Components;
00284     int GlomType;
00287     int StorageType;
00289     int Source;
00291     int Status;
00294     std::vector<vtkStdString> OriginalNames;
00297     std::vector<int> OriginalIndices;
00306     std::vector<int> ObjectTruth;
00308     void Reset();
00309   };
00310 
00312   struct VTKIOEXODUS_EXPORT ObjectInfoType {
00314     int Size;
00316     int Status;
00318     int Id;
00320     vtkStdString Name;
00321   };
00322 
00324   struct VTKIOEXODUS_EXPORT MapInfoType : public ObjectInfoType {
00325   };
00326 
00329   struct VTKIOEXODUS_EXPORT BlockSetInfoType : public ObjectInfoType {
00331     vtkIdType FileOffset;
00336     std::map<vtkIdType,vtkIdType> PointMap;
00341     std::map<vtkIdType,vtkIdType> ReversePointMap;
00345     vtkIdType NextSqueezePoint;
00347     vtkUnstructuredGrid* CachedConnectivity;
00348 
00349     BlockSetInfoType(){this->CachedConnectivity=0;}
00350     BlockSetInfoType(const BlockSetInfoType& block);
00351     ~BlockSetInfoType();
00352     BlockSetInfoType& operator=(const BlockSetInfoType& block);
00353   };
00354 
00356   struct VTKIOEXODUS_EXPORT BlockInfoType : public BlockSetInfoType {
00357     vtkStdString OriginalName; // useful to reset the name if XML metadata is invalid.
00358     vtkStdString TypeName;
00359     // number of boundaries per entry
00360     // The index is the dimensionality of the entry. 0=node, 1=edge, 2=face
00361     int BdsPerEntry[3];
00362     int AttributesPerEntry;
00363     std::vector<vtkStdString> AttributeNames;
00364     std::vector<int> AttributeStatus;
00365     // VTK cell type (a function of TypeName and BdsPerEntry...)
00366     int CellType;
00367     // Number of points per cell as used by VTK
00368     // -- not what's in the file (i.e., BdsPerEntry[0] >= PointsPerCell)
00369     int PointsPerCell;
00370   };
00371 
00373   struct PartInfoType : public ObjectInfoType {
00374     std::vector<int> BlockIndices;
00375   };
00376   struct AssemblyInfoType : public ObjectInfoType {
00377     std::vector<int> BlockIndices;
00378   };
00379   struct MaterialInfoType : public ObjectInfoType {
00380     std::vector<int> BlockIndices;
00381   };
00382 
00384   struct SetInfoType : public BlockSetInfoType {
00385     int DistFact;     // Number of distribution factors
00386                       // (for the entire block, not per array or entry)
00387   };
00388 
00391   enum GlomTypes {
00392     Scalar=0,          
00393     Vector2=1,         
00394     Vector3=2,         
00395     SymmetricTensor=3, 
00396                        //   (order xx, yy, zz, xy, yz, zx)
00397     IntegrationPoint=4 
00398   };
00399 
00401   enum ArraySourceTypes {
00402     Result=0,        
00403                      //   (that vary over time)
00404     Attribute=1,     
00405                      //   (constants over time)
00406     Map=2,           
00407     Generated=3      
00408   };
00409 
00411   vtkTimeStamp InformationTimeStamp;
00412 
00413   friend class vtkExodusIIReader;
00414   friend class vtkPExodusIIReader;
00415 
00416   virtual void SetParser( vtkExodusIIReaderParser* );
00417   vtkGetObjectMacro(Parser,vtkExodusIIReaderParser);
00418 
00419   // Because Parts, Materials, and assemblies are not stored as arrays,
00420   // but rather as maps to the element blocks they make up,
00421   // we cannot use the Get|SetObject__() methods directly.
00422 
00423   int GetNumberOfParts();
00424   const char* GetPartName(int idx);
00425   const char* GetPartBlockInfo(int idx);
00426   int GetPartStatus(int idx);
00427   int GetPartStatus(vtkStdString name);
00428   void SetPartStatus(int idx, int on);
00429   void SetPartStatus(vtkStdString name, int flag);
00430 
00431   int GetNumberOfMaterials();
00432   const char* GetMaterialName(int idx);
00433   int GetMaterialStatus(int idx);
00434   int GetMaterialStatus(vtkStdString name);
00435   void SetMaterialStatus(int idx, int on);
00436   void SetMaterialStatus(vtkStdString name, int flag);
00437 
00438   int GetNumberOfAssemblies();
00439   const char* GetAssemblyName(int idx);
00440   int GetAssemblyStatus(int idx);
00441   int GetAssemblyStatus(vtkStdString name);
00442   void SetAssemblyStatus(int idx, int on);
00443   void SetAssemblyStatus(vtkStdString name, int flag);
00444 
00445   void SetFastPathObjectType(vtkExodusIIReader::ObjectType type)
00446     {this->FastPathObjectType = type;};
00447   void SetFastPathObjectId(vtkIdType id){this->FastPathObjectId = id;};
00448   vtkSetStringMacro(FastPathIdType);
00449 
00450   bool IsXMLMetadataValid();
00451 
00459   void GetInitialObjectStatus( int otype, ObjectInfoType *info );
00460 
00468   void GetInitialObjectArrayStatus( int otype, ArrayInfoType *info );
00469 
00476   void SetInitialObjectStatus( int otype, const char *name, int stat );
00477 
00483   void SetInitialObjectArrayStatus( int otype, const char *name, int stat );
00484 
00485   int UpdateTimeInformation();
00486 
00487   bool ProducedFastPathOutput;
00488 
00489 protected:
00490   vtkExodusIIReaderPrivate();
00491   ~vtkExodusIIReaderPrivate();
00492 
00494   void BuildSIL();
00495 
00498   int VerifyIntegrationPointGlom( int nn,
00499                                   char** np,
00500                                   vtksys::RegularExpression& re,
00501                                   vtkStdString& field,
00502                                   vtkStdString& ele );
00503 
00505   void GlomArrayNames( int i,
00506                        int num_obj,
00507                        int num_vars,
00508                        char** var_names,
00509                        int* truth_tab );
00510 
00512   void PrepareGeneratedArrayInfo();
00513 
00529   int AssembleOutputConnectivity( vtkIdType timeStep,
00530     int otyp, int oidx, int conntypidx, BlockSetInfoType* bsinfop,
00531     vtkUnstructuredGrid* output );
00538   int AssembleOutputPoints( vtkIdType timeStep,
00539     BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output );
00543   int AssembleOutputPointArrays( vtkIdType timeStep,
00544     BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output );
00548   int AssembleOutputCellArrays( vtkIdType timeStep,
00549     int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output );
00554   int AssembleOutputProceduralArrays( vtkIdType timeStep,
00555     int otyp, int oidx, vtkUnstructuredGrid* output );
00557   int AssembleOutputGlobalArrays( vtkIdType timeStep,
00558     int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output );
00565   int AssembleOutputPointMaps( vtkIdType timeStep,
00566     BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output );
00567   int AssembleOutputCellMaps( vtkIdType timeStep,
00568     int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output );
00571   int AssembleArraysOverTime(vtkMultiBlockDataSet* output);
00572 
00574   void InsertBlockCells(
00575     int otyp, int obj, int conn_type, int timeStep, BlockInfoType* binfop );
00576 
00578   void InsertSetCells(
00579     int otyp, int obj, int conn_type, int timeStep, SetInfoType* sinfop );
00580 
00582   void AddPointArray(
00583     vtkDataArray* src, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output );
00584 
00586   void InsertSetNodeCopies(
00587     vtkIntArray* refs, int otyp, int obj, SetInfoType* sinfo );
00588 
00590   void InsertSetCellCopies(
00591     vtkIntArray* refs, int otyp, int obj, SetInfoType* sinfo );
00592 
00594   void InsertSetSides(
00595     vtkIntArray* refs, int otyp, int obj, SetInfoType* sinfo );
00596 
00602   vtkDataArray* GetCacheOrRead( vtkExodusIICacheKey );
00603 
00608   int GetConnTypeIndexFromConnType( int ctyp );
00609 
00614   int GetObjectTypeIndexFromObjectType( int otyp );
00615 
00621   int GetNumberOfObjectsAtTypeIndex( int typeIndex );
00622 
00630   ObjectInfoType* GetObjectInfo( int typeIndex, int objectIndex );
00631 
00638   ObjectInfoType* GetSortedObjectInfo( int objectType, int objectIndex );
00639 
00646   ObjectInfoType* GetUnsortedObjectInfo( int objectType, int objectIndex );
00647 
00652   int GetBlockIndexFromFileGlobalId( int otyp, int refId );
00653 
00658   BlockInfoType* GetBlockFromFileGlobalId( int otyp, int refId );
00659 
00663   vtkIdType GetSqueezePointId( BlockSetInfoType* bsinfop, int i );
00664 
00666   void DetermineVtkCellType( BlockInfoType& binfo );
00667 
00671   ArrayInfoType* FindArrayInfoByName( int otyp, const char* name );
00672 
00676   int IsObjectTypeBlock( int otyp );
00677   int IsObjectTypeSet( int otyp );
00678   int IsObjectTypeMap( int otyp );
00679 
00683   int GetObjectTypeFromMapType( int mtyp );
00684   int GetMapTypeFromObjectType( int otyp );
00685   int GetTemporalTypeFromObjectType( int otyp );
00686 
00690   int GetSetTypeFromSetConnType( int sctyp );
00691 
00695   int GetBlockConnTypeFromBlockType( int btyp );
00696 
00702   void RemoveBeginningAndTrailingSpaces( int len, char **names );
00703 
00705   void ClearConnectivityCaches();
00706 
00710   std::map<int,std::vector<BlockInfoType> > BlockInfo;
00714   std::map<int,std::vector<SetInfoType> > SetInfo;
00720   std::map<int,std::vector<MapInfoType> > MapInfo;
00721 
00722   std::vector<PartInfoType> PartInfo;
00723   std::vector<MaterialInfoType> MaterialInfo;
00724   std::vector<AssemblyInfoType> AssemblyInfo;
00725 
00730   std::map<int,std::vector<int> > SortedObjectIndices;
00732   //  defined on that type.
00733   std::map<int,std::vector<ArrayInfoType> > ArrayInfo;
00734 
00739   std::map<int,std::vector<ArrayInfoType> > InitialArrayInfo;
00740 
00745   std::map<int,std::vector<ObjectInfoType> > InitialObjectInfo;
00746 
00748   int AppWordSize;
00749   int DiskWordSize;
00750 
00754   float ExodusVersion;
00755 
00757   int Exoid;
00758 
00760   struct ex_init_params ModelParameters;
00761 
00763   std::vector<double> Times;
00764 
00766   int TimeStep;
00767 
00771   double ModeShapeTime;
00772 
00773   int GenerateObjectIdArray;
00774   int GenerateGlobalIdArray;
00775   int GenerateFileIdArray;
00776   int GenerateGlobalElementIdArray;
00777   int GenerateGlobalNodeIdArray;
00778   int GenerateImplicitElementIdArray;
00779   int GenerateImplicitNodeIdArray;
00780 
00784   int FileId;
00785 
00787   vtkExodusIICache* Cache;
00788   //
00790   double CacheSize;
00791 
00792   int ApplyDisplacements;
00793   float DisplacementMagnitude;
00794   int HasModeShapes;
00795   int AnimateModeShapes;
00796 
00808   int SqueezePoints;
00809 
00813   vtkExodusIIReader* Parent;
00814 
00815   vtkExodusIIReaderParser* Parser;
00816 
00817   vtkExodusIIReader::ObjectType FastPathObjectType;
00818   vtkIdType FastPathObjectId;
00819   char* FastPathIdType;
00820 
00821   vtkMutableDirectedGraph* SIL;
00822 private:
00823   vtkExodusIIReaderPrivate( const vtkExodusIIReaderPrivate& ); // Not implemented.
00824   void operator = ( const vtkExodusIIReaderPrivate& ); // Not implemented.
00825 };
00826 
00827 #endif // __vtkExodusIIReaderPrivate_h
00828 // VTK-HeaderTest-Exclude: vtkExodusIIReaderPrivate.h