VTK  9.0.1
vtkExodusIIReaderPrivate.h
Go to the documentation of this file.
1 #ifndef vtkExodusIIReaderPrivate_h
2 #define vtkExodusIIReaderPrivate_h
3 #ifndef __VTK_WRAP__
4 #ifndef VTK_WRAPPING_CXX
5 
6 // Do not include this file directly. It is only for use
7 // from inside the ExodusII reader and its descendants.
8 
9 #include "vtkExodusIICache.h"
10 #include "vtkExodusIIReader.h"
11 #include "vtkStdString.h"
12 #include "vtkToolkits.h" // make sure VTK_USE_PARALLEL is properly set
13 #include "vtksys/RegularExpression.hxx"
14 
15 #include <map>
16 #include <vector>
17 
18 #include "vtkIOExodusModule.h" // For export macro
19 #include "vtk_exodusII.h"
21 class vtkIdTypeArray;
24 class vtkTypeInt64Array;
26 
30 class VTKIOEXODUS_EXPORT vtkExodusIIReaderPrivate : public vtkObject
31 {
32 public:
33  static vtkExodusIIReaderPrivate* New();
34  void PrintData(ostream& os, vtkIndent indent);
36  // virtual void Modified();
37 
39  int OpenFile(const char* filename);
40 
42  int CloseFile();
43 
45  int RequestInformation();
46 
48  vtkMutableDirectedGraph* GetSIL() { return this->SIL; }
49 
51  int RequestData(vtkIdType timeStep, vtkMultiBlockDataSet* output);
52 
57  int SetUpEmptyGrid(vtkMultiBlockDataSet* output);
58 
70  void Reset();
71 
76  void ResetSettings();
77 
79  void ResetCache();
80 
82  void SetCacheSize(double size);
83 
85  vtkGetMacro(CacheSize, double);
86 
91  int GetNumberOfTimeSteps() { return (int)this->Times.size(); }
92 
95  vtkGetMacro(SqueezePoints, int);
96 
99  void SetSqueezePoints(int sp);
100 
103  vtkBooleanMacro(SqueezePoints, int);
104 
106  int GetNumberOfNodes();
107 
112  int GetNumberOfObjectsOfType(int otype);
113 
124  int GetNumberOfObjectArraysOfType(int otype);
125 
130  const char* GetObjectName(int otype, int i);
131 
136  int GetObjectId(int otype, int i);
137 
144  int GetObjectSize(int otype, int i);
145 
150  int GetObjectStatus(int otype, int i);
151 
157  int GetUnsortedObjectStatus(int otype, int i);
158 
163  void SetObjectStatus(int otype, int i, int stat);
164 
170  void SetUnsortedObjectStatus(int otype, int i, int stat);
171 
176  const char* GetObjectArrayName(int otype, int i);
177 
182  int GetNumberOfObjectArrayComponents(int otype, int i);
183 
188  int GetObjectArrayStatus(int otype, int i);
189 
194  void SetObjectArrayStatus(int otype, int i, int stat);
195 
202  int GetNumberOfObjectAttributes(int objectType, int objectIndex);
203  const char* GetObjectAttributeName(int objectType, int objectIndex, int attributeIndex);
204  int GetObjectAttributeIndex(int objectType, int objectIndex, const char* attribName);
205  int GetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex);
206  void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status);
207 
209  vtkGetMacro(GenerateObjectIdArray, vtkTypeBool);
210  vtkSetMacro(GenerateObjectIdArray, vtkTypeBool);
211  static const char* GetObjectIdArrayName() { return "ObjectId"; }
212 
213  vtkSetMacro(GenerateGlobalElementIdArray, vtkTypeBool);
214  vtkGetMacro(GenerateGlobalElementIdArray, vtkTypeBool);
215  static const char* GetGlobalElementIdArrayName() { return "GlobalElementId"; }
216 
217  vtkSetMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
218  vtkGetMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
219  static const char* GetGlobalNodeIdArrayName() { return "GlobalNodeId"; }
220 
221  vtkSetMacro(GenerateImplicitElementIdArray, vtkTypeBool);
222  vtkGetMacro(GenerateImplicitElementIdArray, vtkTypeBool);
223  static const char* GetImplicitElementIdArrayName() { return "ImplicitElementId"; }
224 
225  vtkSetMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
226  vtkGetMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
227  static const char* GetImplicitNodeIdArrayName() { return "ImplicitNodeId"; }
228 
232  vtkSetMacro(GenerateFileIdArray, vtkTypeBool);
233  vtkGetMacro(GenerateFileIdArray, vtkTypeBool);
234  static const char* GetFileIdArrayName() { return "FileId"; }
235 
237  vtkSetMacro(FileId, int);
238  vtkGetMacro(FileId, int);
239 
240  static const char* GetGlobalVariableValuesArrayName() { return "GlobalVariableValues"; }
241  static const char* GetGlobalVariableNamesArrayName() { return "GlobalVariableNames"; }
242 
243  virtual void SetApplyDisplacements(vtkTypeBool d);
244  vtkGetMacro(ApplyDisplacements, vtkTypeBool);
245 
246  virtual void SetDisplacementMagnitude(double s);
247  vtkGetMacro(DisplacementMagnitude, double);
248 
249  vtkSetMacro(HasModeShapes, int);
250  vtkGetMacro(HasModeShapes, int);
251 
252  vtkSetMacro(ModeShapeTime, double);
253  vtkGetMacro(ModeShapeTime, double);
254 
255  vtkSetMacro(AnimateModeShapes, int);
256  vtkGetMacro(AnimateModeShapes, int);
257 
258  vtkSetMacro(IgnoreFileTime, bool);
259  vtkGetMacro(IgnoreFileTime, bool);
260 
261  vtkDataArray* FindDisplacementVectors(int timeStep);
262 
263  const struct ex_init_params* GetModelParams() const { return &this->ModelParameters; }
264 
266  struct VTKIOEXODUS_EXPORT ArrayInfoType
267  {
278  int GlomType;
283  int Source;
285  int Status;
288  std::vector<vtkStdString> OriginalNames;
291  std::vector<int> OriginalIndices;
300  std::vector<int> ObjectTruth;
302  void Reset();
303  };
304 
306  struct VTKIOEXODUS_EXPORT ObjectInfoType
307  {
309  int Size;
311  int Status;
313  int Id;
316  };
317 
319  struct VTKIOEXODUS_EXPORT MapInfoType : public ObjectInfoType
320  {
321  };
322 
325  struct VTKIOEXODUS_EXPORT BlockSetInfoType : public ObjectInfoType
326  {
333  std::map<vtkIdType, vtkIdType> PointMap;
338  std::map<vtkIdType, vtkIdType> ReversePointMap;
345 
346  BlockSetInfoType() { this->CachedConnectivity = nullptr; }
347  BlockSetInfoType(const BlockSetInfoType& block);
348  ~BlockSetInfoType();
349  BlockSetInfoType& operator=(const BlockSetInfoType& block);
350  };
351 
353  struct VTKIOEXODUS_EXPORT BlockInfoType : public BlockSetInfoType
354  {
355  vtkStdString OriginalName; // useful to reset the name if XML metadata is invalid.
357  // number of boundaries per entry
358  // The index is the dimensionality of the entry. 0=node, 1=edge, 2=face
359  int BdsPerEntry[3];
361  std::vector<vtkStdString> AttributeNames;
362  std::vector<int> AttributeStatus;
363  // VTK cell type (a function of TypeName and BdsPerEntry...)
364  int CellType;
365  // Number of points per cell as used by VTK
366  // -- not what's in the file (i.e., BdsPerEntry[0] >= PointsPerCell)
368  };
369 
371  struct PartInfoType : public ObjectInfoType
372  {
373  std::vector<int> BlockIndices;
374  };
376  {
377  std::vector<int> BlockIndices;
378  };
380  {
381  std::vector<int> BlockIndices;
382  };
383 
386  {
387  int DistFact; // Number of distribution factors
388  // (for the entire block, not per array or entry)
389  };
390 
394  {
395  Scalar = 0,
396  Vector2 = 1,
397  Vector3 = 2,
398  SymmetricTensor = 3,
399  // (order xx, yy, zz, xy, yz, zx)
400  IntegrationPoint = 4
401  };
402 
405  {
406  Result = 0,
407  // (that vary over time)
408  Attribute = 1,
409  // (constants over time)
410  Map = 2,
411  Generated = 3
412  };
413 
416 
417  friend class vtkExodusIIReader;
418  friend class vtkPExodusIIReader;
419 
420  virtual void SetParser(vtkExodusIIReaderParser*);
421  vtkGetObjectMacro(Parser, vtkExodusIIReaderParser);
422 
423  // BUG #15632: This method allows vtkPExodusIIReader to pass time information
424  // from one spatial file to another and avoiding have to read it for each of
425  // the files.
426  void SetTimesOverrides(const std::vector<double>& times)
427  {
428  this->Times = times;
429  this->SkipUpdateTimeInformation = true;
430  }
431 
432  // Because Parts, Materials, and assemblies are not stored as arrays,
433  // but rather as maps to the element blocks they make up,
434  // we cannot use the Get|SetObject__() methods directly.
435 
436  int GetNumberOfParts();
437  const char* GetPartName(int idx);
438  const char* GetPartBlockInfo(int idx);
439  int GetPartStatus(int idx);
440  int GetPartStatus(const vtkStdString& name);
441  void SetPartStatus(int idx, int on);
442  void SetPartStatus(const vtkStdString& name, int flag);
443 
444  int GetNumberOfMaterials();
445  const char* GetMaterialName(int idx);
446  int GetMaterialStatus(int idx);
447  int GetMaterialStatus(const vtkStdString& name);
448  void SetMaterialStatus(int idx, int on);
449  void SetMaterialStatus(const vtkStdString& name, int flag);
450 
451  int GetNumberOfAssemblies();
452  const char* GetAssemblyName(int idx);
453  int GetAssemblyStatus(int idx);
454  int GetAssemblyStatus(const vtkStdString& name);
455  void SetAssemblyStatus(int idx, int on);
456  void SetAssemblyStatus(const vtkStdString& name, int flag);
457 
459  {
460  this->FastPathObjectType = type;
461  }
462  void SetFastPathObjectId(vtkIdType id) { this->FastPathObjectId = id; }
463  vtkSetStringMacro(FastPathIdType);
464 
465  bool IsXMLMetadataValid();
466 
474  void GetInitialObjectStatus(int otype, ObjectInfoType* info);
475 
483  void GetInitialObjectArrayStatus(int otype, ArrayInfoType* info);
484 
491  void SetInitialObjectStatus(int otype, const char* name, int stat);
492 
498  void SetInitialObjectArrayStatus(int otype, const char* name, int stat);
499 
500  int UpdateTimeInformation();
501 
503 
504 protected:
506  ~vtkExodusIIReaderPrivate() override;
507 
509  void BuildSIL();
510 
514  int nn, char** np, vtksys::RegularExpression& re, vtkStdString& field, vtkStdString& ele);
515 
517  void GlomArrayNames(int i, int num_obj, int num_vars, char** var_names, int* truth_tab);
518 
521 
538  int AssembleOutputConnectivity(vtkIdType timeStep, int otyp, int oidx, int conntypidx,
539  BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
547  vtkIdType timeStep, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
552  vtkIdType timeStep, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
557  vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
563  vtkIdType timeStep, int otyp, int oidx, vtkUnstructuredGrid* output);
566  vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
574  vtkIdType timeStep, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
576  vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
580 
597  vtkIdType GetPolyhedronFaceConnectivity(vtkIdType fileLocalFaceId, vtkIdType*& facePtIds);
598 
601 
604  BlockInfoType* binfo, vtkIntArray* facesPerCell, vtkIdTypeArray* exoCellConn);
605 
607  void InsertBlockCells(int otyp, int obj, int conn_type, int timeStep, BlockInfoType* binfop);
608 
610  void InsertSetCells(int otyp, int obj, int conn_type, int timeStep, SetInfoType* sinfop);
611 
613  void AddPointArray(vtkDataArray* src, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
614 
616  void InsertSetNodeCopies(vtkIdTypeArray* refs, int otyp, int obj, SetInfoType* sinfo);
617 
619  void InsertSetCellCopies(vtkIdTypeArray* refs, int otyp, int obj, SetInfoType* sinfo);
620 
622  void InsertSetSides(vtkIdTypeArray* refs, int otyp, int obj, SetInfoType* sinfo);
623 
630 
635  int GetConnTypeIndexFromConnType(int ctyp);
636 
641  int GetObjectTypeIndexFromObjectType(int otyp);
642 
648  int GetNumberOfObjectsAtTypeIndex(int typeIndex);
649 
657  ObjectInfoType* GetObjectInfo(int typeIndex, int objectIndex);
658 
665  ObjectInfoType* GetSortedObjectInfo(int objectType, int objectIndex);
666 
673  ObjectInfoType* GetUnsortedObjectInfo(int objectType, int objectIndex);
674 
679  int GetBlockIndexFromFileGlobalId(int otyp, int refId);
680 
685  BlockInfoType* GetBlockFromFileGlobalId(int otyp, int refId);
686 
691 
693  void DetermineVtkCellType(BlockInfoType& binfo);
694 
698  ArrayInfoType* FindArrayInfoByName(int otyp, const char* name);
699 
703  int IsObjectTypeBlock(int otyp);
704  int IsObjectTypeSet(int otyp);
705  int IsObjectTypeMap(int otyp);
706 
710  int GetObjectTypeFromMapType(int mtyp);
711  int GetMapTypeFromObjectType(int otyp);
712  int GetTemporalTypeFromObjectType(int otyp);
713 
717  int GetSetTypeFromSetConnType(int sctyp);
718 
722  int GetBlockConnTypeFromBlockType(int btyp);
723 
729  void RemoveBeginningAndTrailingSpaces(int len, char** names, int maxNameLength);
730 
733 
737  std::map<int, std::vector<BlockInfoType> > BlockInfo;
741  std::map<int, std::vector<SetInfoType> > SetInfo;
747  std::map<int, std::vector<MapInfoType> > MapInfo;
748 
749  std::vector<PartInfoType> PartInfo;
750  std::vector<MaterialInfoType> MaterialInfo;
751  std::vector<AssemblyInfoType> AssemblyInfo;
752 
757  std::map<int, std::vector<int> > SortedObjectIndices;
759  // defined on that type.
760  std::map<int, std::vector<ArrayInfoType> > ArrayInfo;
761 
766  std::map<int, std::vector<ArrayInfoType> > InitialArrayInfo;
767 
772  std::map<int, std::vector<ObjectInfoType> > InitialObjectInfo;
773 
777 
782 
784  int Exoid;
785 
787  struct ex_init_params ModelParameters;
788 
790  std::vector<double> Times;
792 
797 
805 
809  int FileId;
810 
813  //
815  double CacheSize;
816 
821 
823 
836 
841 
843 
852  std::map<int, std::vector<std::vector<vtkIdType> > > PolyhedralFaceConnArrays;
853 
857 
859 
860 private:
862  void operator=(const vtkExodusIIReaderPrivate&) = delete;
863 };
864 
865 #endif
866 #endif
867 #endif // vtkExodusIIReaderPrivate_h
868 // VTK-HeaderTest-Exclude: vtkExodusIIReaderPrivate.h
void InsertSetCells(int otyp, int obj, int conn_type, int timeStep, SetInfoType *sinfop)
Insert cells from a specified set into a mesh.
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:34
std::vector< int > OriginalIndices
The index of each component of the array as ordered by the Exodus file.
std::vector< MaterialInfoType > MaterialInfo
int Components
The number of components in the array.
static const char * GetObjectIdArrayName()
vtkMutableDirectedGraph * SIL
const char * GetAssemblyName(int idx)
void SetFastPathObjectType(vtkExodusIIReader::ObjectType type)
void RemoveBeginningAndTrailingSpaces(int len, char **names, int maxNameLength)
Function to trim space from names retrieved with ex_get_var_names.
const char * GetPartName(int idx)
abstract base class for most VTK objects
Definition: vtkObject.h:62
A struct to hold information about Exodus objects (blocks, sets, maps)
std::map< int, std::vector< ObjectInfoType > > InitialObjectInfo
Maps an object type (EX_ELEM_BLOCK, EX_NODE_SET, ...) to a list of objects defined on that type...
std::map< int, std::vector< ArrayInfoType > > ArrayInfo
Maps an object type (EX_ELEM_BLOCK, EX_NODE_SET, ...) to a list of arrays.
void SetFastPathObjectId(vtkIdType id)
int GetNumberOfTimeSteps()
Return the number of time steps in the open file.
void GetInitialObjectStatus(int otype, ObjectInfoType *info)
For a given object type, looks for an object in the collection of initial objects of the same name...
void InsertSetNodeCopies(vtkIdTypeArray *refs, int otyp, int obj, SetInfoType *sinfo)
Insert cells referenced by a node set.
void PrepareGeneratedArrayInfo()
Add generated array information to array info lists.
int GetMapTypeFromObjectType(int otyp)
int AssembleOutputPoints(vtkIdType timeStep, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Fill the output grid&#39;s point coordinates array.
int GetObjectTypeFromMapType(int mtyp)
Given a map type (NODE_MAP, EDGE_MAP, ...) return the associated object type (NODAL, EDGE_BLOCK, ...) or vice-versa.
int GetNumberOfObjectsAtTypeIndex(int typeIndex)
Return the number of objects of the given type.
record modification and/or execution time
Definition: vtkTimeStamp.h:32
std::vector< int > ObjectTruth
A map describing which objects the variable is defined on.
A struct to hold information about Exodus maps.
void SetTimesOverrides(const std::vector< double > &times)
std::map< int, std::vector< int > > SortedObjectIndices
Maps an object type to vector of indices that reorder objects of that type by their IDs...
void FreePolyhedronFaceArrays()
Free any arrays held by PolyhedralFaceConnArrays (for polyhedral-face-connectivity lookup)...
BlockInfoType * GetBlockFromFileGlobalId(int otyp, int refId)
Get the block containing the entity referenced by the specified file-global ID.
vtkExodusIIReader * Parent
Pointer to owning reader...
void DetermineVtkCellType(BlockInfoType &binfo)
Determine the VTK cell type for a given edge/face/element block.
int AssembleOutputPointArrays(vtkIdType timeStep, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Add the requested arrays to the output grid&#39;s point data.
dynamic, self-adjusting array of vtkIdType
void ClearConnectivityCaches()
Delete any cached connectivity information (for all blocks and sets)
ObjectType
Extra cell data array that can be generated.
std::map< int, std::vector< std::vector< vtkIdType > > > PolyhedralFaceConnArrays
Face connectivity for polyhedra.
int vtkIdType
Definition: vtkType.h:338
This class holds metadata for an Exodus file.
vtkUnstructuredGrid * CachedConnectivity
Cached cell connectivity arrays for mesh.
vtkMutableDirectedGraph * GetSIL()
Returns the SIL. This valid only after BuildSIL() has been called.
void GlomArrayNames(int i, int num_obj, int num_vars, char **var_names, int *truth_tab)
Aggregate Exodus array names into VTK arrays with multiple components.
int GetBlockConnTypeFromBlockType(int btyp)
Given a block type (EDGE_BLOCK, ...), return the associated block connectivity type (EDGE_BLOCK_CONN...
void InsertSetSides(vtkIdTypeArray *refs, int otyp, int obj, SetInfoType *sinfo)
Insert cells referenced by a side set.
static const char * GetGlobalVariableValuesArrayName()
int vtkTypeBool
Definition: vtkABI.h:69
double CacheSize
The size of the cache in MiB.
std::map< int, std::vector< ArrayInfoType > > InitialArrayInfo
Maps an object type (EX_ELEM_BLOCK, EX_NODE_SET, ...) to a list of arrays defined on that type...
void InsertBlockPolyhedra(BlockInfoType *binfo, vtkIntArray *facesPerCell, vtkIdTypeArray *exoCellConn)
Insert polyhedral cells (called from InsertBlockCells when a block is polyhedral).
vtkTimeStamp InformationTimeStamp
Time stamp from last time we were in RequestInformation.
void GetInitialObjectArrayStatus(int otype, ArrayInfoType *info)
For a given array type, looks for an object in the collection of initial objects of the same name...
int Source
The source of the array (Result or Attribute)
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:39
int AssembleOutputProceduralArrays(vtkIdType timeStep, int otyp, int oidx, vtkUnstructuredGrid *output)
Add procedurally generated arrays to an output mesh.
ObjectInfoType * GetObjectInfo(int typeIndex, int objectIndex)
Return a pointer to the ObjectInfo of the specified type and index.
int AssembleOutputGlobalArrays(vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Add mesh-global field data such as QA records to the output mesh.
void BuildSIL()
Build SIL. This must be called only after RequestInformation().
static const char * GetFileIdArrayName()
int StorageType
Storage type of array (a type that can be passed to vtkDataArray::Create())
double ModeShapeTime
The time value.
std::vector< AssemblyInfoType > AssemblyInfo
vtkStdString Name
The name of the array.
a simple class to control print indentation
Definition: vtkIndent.h:33
float ExodusVersion
The version of Exodus that wrote the currently open file (or a negative number otherwise).
std::map< int, std::vector< MapInfoType > > MapInfo
Maps a map type (EX_ELEM_MAP, ..., EX_NODE_MAP) to a list of maps of that type.
Read Exodus II files (.exii)
int GetSetTypeFromSetConnType(int sctyp)
Given a set connectivity type (NODE_SET_CONN, ...), return the associated object type (NODE_SET...
GlomTypes
Tags to indicate how single-component Exodus arrays are glommed (aggregated) into multi-component VTK...
dataset represents arbitrary combinations of all possible cell types
vtkIdType FileOffset
Id (1-based) of first entry in file-local list across all blocks in file.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:49
static const char * GetGlobalElementIdArrayName()
void SetAssemblyStatus(int idx, int on)
std::map< vtkIdType, vtkIdType > PointMap
A map from nodal IDs in an Exodus file to nodal IDs in the output mesh.
vtkIdType GetPolyhedronFaceConnectivity(vtkIdType fileLocalFaceId, vtkIdType *&facePtIds)
Fetch the face-connectivity for one face of one polyhedron.
An editable directed graph.
int IsObjectTypeBlock(int otyp)
Does the specified object type match? Avoid using these...
std::map< int, std::vector< SetInfoType > > SetInfo
Maps a set type (EX_ELEM_SET, ..., EX_NODE_SET) to a list of sets of that type.
int GetPartStatus(int idx)
vtkExodusIIReader::ObjectType FastPathObjectType
int Status
Whether or not the array should be loaded by RequestData.
int GetMaterialStatus(int idx)
vtkIdType GetSqueezePointId(BlockSetInfoType *bsinfop, int i)
Find or create a new SqueezePoint ID (unique sequential list of points referenced by cells in blocks/...
virtual void SetParser(vtkExodusIIReaderParser *)
~vtkExodusIIReaderPrivate() override
int IsObjectTypeSet(int otyp)
A struct to hold information about Exodus blocks.
int VerifyIntegrationPointGlom(int nn, char **np, vtksys::RegularExpression &re, vtkStdString &field, vtkStdString &ele)
Returns true when order and text of names are consistent with integration points. ...
int Id
User-assigned identification number.
int AssembleOutputPointMaps(vtkIdType timeStep, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Add maps to an output mesh.
int GetAssemblyStatus(int idx)
A struct to hold information about Exodus blocks or sets (they have some members in common) ...
static const char * GetImplicitElementIdArrayName()
std::map< int, std::vector< BlockInfoType > > BlockInfo
Maps a block type (EX_ELEM_BLOCK, EX_FACE_BLOCK, ...) to a list of blocks of that type...
internal parser used by vtkExodusIIReader.
int AssembleOutputCellArrays(vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Add the requested arrays to the output grid&#39;s cell data.
Read exodus 2 files .ex2.
std::vector< vtkStdString > OriginalNames
The name of each component of the array as defined by the Exodus file.
vtkExodusIICache * Cache
A least-recently-used cache to hold raw arrays.
void SetInitialObjectArrayStatus(int otype, const char *name, int stat)
For a given array type, creates and stores an ArrayInfoType object using the given name and status...
std::map< vtkIdType, vtkIdType > ReversePointMap
A map from nodal ids in the output mesh to those in an Exodus file.
ObjectInfoType * GetUnsortedObjectInfo(int objectType, int objectIndex)
Return a pointer to the ObjectInfo of the specified type and index, but using indices sorted by objec...
int GetTemporalTypeFromObjectType(int otyp)
int GetConnTypeIndexFromConnType(int ctyp)
Return the index of an object type (in a private list of all object types).
Composite dataset that organizes datasets into blocks.
static const char * GetGlobalNodeIdArrayName()
void SetPartStatus(int idx, int on)
const struct ex_init_params * GetModelParams() const
int Status
Should the reader load this block?
A struct to hold information about Exodus sets.
static const char * GetImplicitNodeIdArrayName()
ObjectInfoType * GetSortedObjectInfo(int objectType, int objectIndex)
Return a pointer to the ObjectInfo of the specified type and index, but using indices sorted by objec...
A struct to hold information about Exodus blocks.
vtkDataArray * GetCacheOrRead(vtkExodusIICacheKey)
Return an array for the specified cache key.
int GetObjectTypeIndexFromObjectType(int otyp)
Return the index of an object type (in a private list of all object types).
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
int AppWordSize
These aren&#39;t the variables you&#39;re looking for.
int AssembleOutputCellMaps(vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
void AddPointArray(vtkDataArray *src, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Add a point array to an output grid&#39;s point data, squeezing if necessary.
int SqueezePoints
Should the reader output only points used by elements in the output mesh, or all the points...
ArraySourceTypes
Tags to indicate the source of values for an array.
int Exoid
The handle of the currently open file.
const char * GetPartBlockInfo(int idx)
std::vector< double > Times
A list of time steps for which results variables are stored.
vtkIdType NextSqueezePoint
The next vtk ID to use for a connectivity entry when point squeezing is on and no point ID exists...
int Size
Number of entries in this block.
int IsObjectTypeMap(int otyp)
std::vector< PartInfoType > PartInfo
vtkExodusIIReaderParser * Parser
ArrayInfoType * FindArrayInfoByName(int otyp, const char *name)
Find an ArrayInfo object for a specific object type using the name as a key.
void InsertBlockCells(int otyp, int obj, int conn_type, int timeStep, BlockInfoType *binfop)
Insert cells from a specified block into a mesh.
void SetMaterialStatus(int idx, int on)
A struct to hold information about time-varying arrays.
int AssembleOutputConnectivity(vtkIdType timeStep, int otyp, int oidx, int conntypidx, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Read connectivity information and populate an unstructured grid with cells corresponding to a single ...
int GetBlockIndexFromFileGlobalId(int otyp, int refId)
Get the index of the block containing the entity referenced by the specified file-global ID...
BlockSetInfoType & operator=(const BlockSetInfoType &block)
int GlomType
The type of "glomming" performed.
static const char * GetGlobalVariableNamesArrayName()
void SetInitialObjectStatus(int otype, const char *name, int stat)
For a given object type, creates and stores an ObjectInfoType object using the given name and status...
void InsertSetCellCopies(vtkIdTypeArray *refs, int otyp, int obj, SetInfoType *sinfo)
Insert cells referenced by an edge, face, or element set.
int AssembleArraysOverTime(vtkMultiBlockDataSet *output)
Add fast-path time-varying data to field data of an output block or set.
const char * GetMaterialName(int idx)