VTK  9.1.0
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 "vtkDeprecation.h" // for deprecation macros
10 #include "vtkExodusIICache.h" // for vtkExodusIICacheKey
11 #include "vtkExodusIIReader.h" // for vtkExodusIIReader
12 #include "vtkObject.h"
13 #include "vtkStdString.h" // for vtkStdString
14 #include "vtksys/RegularExpression.hxx" // for vtksys::RegularExpression
15 
16 #include <map> // for std::map
17 #include <vector> // for std::vector
18 
19 #include "vtkIOExodusModule.h" // For export macro
20 #include "vtk_exodusII.h" // for exodus APIs
21 class vtkDataArray;
23 class vtkIdTypeArray;
26 class vtkTypeInt64Array;
28 
32 class VTKIOEXODUS_EXPORT vtkExodusIIReaderPrivate : public vtkObject
33 {
34 public:
36  VTK_DEPRECATED_IN_9_1_0("Renamed to PrintSelf")
37  void PrintData(ostream& os, vtkIndent indent) { this->PrintSelf(os, indent); }
38  void PrintSelf(ostream& os, vtkIndent indent) override;
40  // virtual void Modified();
41 
43  int OpenFile(const char* filename);
44 
46  int CloseFile();
47 
50 
52  vtkMutableDirectedGraph* GetSIL() { return this->SIL; }
53 
55  int RequestData(vtkIdType timeStep, vtkMultiBlockDataSet* output);
56 
62 
74  void Reset();
75 
80  void ResetSettings();
81 
83  void ResetCache();
84 
86  void SetCacheSize(double size);
87 
89  vtkGetMacro(CacheSize, double);
90 
95  int GetNumberOfTimeSteps() { return (int)this->Times.size(); }
96 
99  vtkGetMacro(SqueezePoints, int);
100 
103  void SetSqueezePoints(int sp);
104 
107  vtkBooleanMacro(SqueezePoints, int);
108 
111 
116  int GetNumberOfObjectsOfType(int otype);
117 
129 
134  const char* GetObjectName(int otype, int i);
135 
140  int GetObjectId(int otype, int i);
141 
148  int GetObjectSize(int otype, int i);
149 
154  int GetObjectStatus(int otype, int i);
155 
161  int GetUnsortedObjectStatus(int otype, int i);
162 
167  void SetObjectStatus(int otype, int i, int stat);
168 
174  void SetUnsortedObjectStatus(int otype, int i, int stat);
175 
180  const char* GetObjectArrayName(int otype, int i);
181 
186  int GetNumberOfObjectArrayComponents(int otype, int i);
187 
192  int GetObjectArrayStatus(int otype, int i);
193 
198  void SetObjectArrayStatus(int otype, int i, int stat);
199 
206  int GetNumberOfObjectAttributes(int objectType, int objectIndex);
207  const char* GetObjectAttributeName(int objectType, int objectIndex, int attributeIndex);
208  int GetObjectAttributeIndex(int objectType, int objectIndex, const char* attribName);
209  int GetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex);
210  void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status);
211 
213  vtkGetMacro(GenerateObjectIdArray, vtkTypeBool);
214  vtkSetMacro(GenerateObjectIdArray, vtkTypeBool);
215  static const char* GetObjectIdArrayName() { return "ObjectId"; }
216 
217  vtkSetMacro(GenerateGlobalElementIdArray, vtkTypeBool);
218  vtkGetMacro(GenerateGlobalElementIdArray, vtkTypeBool);
219  static const char* GetGlobalElementIdArrayName() { return "GlobalElementId"; }
220 
221  vtkSetMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
222  vtkGetMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
223  static const char* GetGlobalNodeIdArrayName() { return "GlobalNodeId"; }
224 
225  vtkSetMacro(GenerateImplicitElementIdArray, vtkTypeBool);
226  vtkGetMacro(GenerateImplicitElementIdArray, vtkTypeBool);
227  static const char* GetImplicitElementIdArrayName() { return "ImplicitElementId"; }
228 
229  vtkSetMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
230  vtkGetMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
231  static const char* GetImplicitNodeIdArrayName() { return "ImplicitNodeId"; }
232 
236  vtkSetMacro(GenerateFileIdArray, vtkTypeBool);
237  vtkGetMacro(GenerateFileIdArray, vtkTypeBool);
238  static const char* GetFileIdArrayName() { return "FileId"; }
239 
241  vtkSetMacro(FileId, int);
242  vtkGetMacro(FileId, int);
243 
244  static const char* GetGlobalVariableValuesArrayName() { return "GlobalVariableValues"; }
245  static const char* GetGlobalVariableNamesArrayName() { return "GlobalVariableNames"; }
246 
248  vtkGetMacro(ApplyDisplacements, vtkTypeBool);
249 
250  virtual void SetDisplacementMagnitude(double s);
251  vtkGetMacro(DisplacementMagnitude, double);
252 
253  vtkSetMacro(HasModeShapes, int);
254  vtkGetMacro(HasModeShapes, int);
255 
256  vtkSetMacro(ModeShapeTime, double);
257  vtkGetMacro(ModeShapeTime, double);
258 
259  vtkSetMacro(AnimateModeShapes, int);
260  vtkGetMacro(AnimateModeShapes, int);
261 
262  vtkSetMacro(IgnoreFileTime, bool);
263  vtkGetMacro(IgnoreFileTime, bool);
264 
266 
267  const struct ex_init_params* GetModelParams() const { return &this->ModelParameters; }
268 
270  struct VTKIOEXODUS_EXPORT ArrayInfoType
271  {
282  int GlomType;
287  int Source;
289  int Status;
292  std::vector<vtkStdString> OriginalNames;
295  std::vector<int> OriginalIndices;
304  std::vector<int> ObjectTruth;
306  void Reset();
307  };
308 
310  struct VTKIOEXODUS_EXPORT ObjectInfoType
311  {
313  int Size;
315  int Status;
317  int Id;
320  };
321 
323  struct VTKIOEXODUS_EXPORT MapInfoType : public ObjectInfoType
324  {
325  };
326 
329  struct VTKIOEXODUS_EXPORT BlockSetInfoType : public ObjectInfoType
330  {
337  std::map<vtkIdType, vtkIdType> PointMap;
342  std::map<vtkIdType, vtkIdType> ReversePointMap;
349 
350  BlockSetInfoType() { this->CachedConnectivity = nullptr; }
354  };
355 
357  struct VTKIOEXODUS_EXPORT BlockInfoType : public BlockSetInfoType
358  {
359  vtkStdString OriginalName; // useful to reset the name if XML metadata is invalid.
361  // number of boundaries per entry
362  // The index is the dimensionality of the entry. 0=node, 1=edge, 2=face
363  int BdsPerEntry[3];
365  std::vector<vtkStdString> AttributeNames;
366  std::vector<int> AttributeStatus;
367  // VTK cell type (a function of TypeName and BdsPerEntry...)
368  int CellType;
369  // Number of points per cell as used by VTK
370  // -- not what's in the file (i.e., BdsPerEntry[0] >= PointsPerCell)
372  };
373 
375  struct PartInfoType : public ObjectInfoType
376  {
377  std::vector<int> BlockIndices;
378  };
380  {
381  std::vector<int> BlockIndices;
382  };
384  {
385  std::vector<int> BlockIndices;
386  };
387 
390  {
391  int DistFact; // Number of distribution factors
392  // (for the entire block, not per array or entry)
393  };
394 
398  {
399  Scalar = 0,
400  Vector2 = 1,
401  Vector3 = 2,
402  SymmetricTensor = 3,
403  // (order xx, yy, zz, xy, yz, zx)
404  IntegrationPoint = 4
405  };
406 
409  {
410  Result = 0,
411  // (that vary over time)
412  Attribute = 1,
413  // (constants over time)
414  Map = 2,
415  Generated = 3
416  };
417 
420 
421  friend class vtkExodusIIReader;
422  friend class vtkPExodusIIReader;
423 
425  vtkGetObjectMacro(Parser, vtkExodusIIReaderParser);
426 
427  // BUG #15632: This method allows vtkPExodusIIReader to pass time information
428  // from one spatial file to another and avoiding have to read it for each of
429  // the files.
430  void SetTimesOverrides(const std::vector<double>& times)
431  {
432  this->Times = times;
433  this->SkipUpdateTimeInformation = true;
434  }
435 
436  // Because Parts, Materials, and assemblies are not stored as arrays,
437  // but rather as maps to the element blocks they make up,
438  // we cannot use the Get|SetObject__() methods directly.
439 
441  const char* GetPartName(int idx);
442  const char* GetPartBlockInfo(int idx);
443  int GetPartStatus(int idx);
445  void SetPartStatus(int idx, int on);
446  void SetPartStatus(const vtkStdString& name, int flag);
447 
449  const char* GetMaterialName(int idx);
450  int GetMaterialStatus(int idx);
452  void SetMaterialStatus(int idx, int on);
453  void SetMaterialStatus(const vtkStdString& name, int flag);
454 
456  const char* GetAssemblyName(int idx);
457  int GetAssemblyStatus(int idx);
459  void SetAssemblyStatus(int idx, int on);
460  void SetAssemblyStatus(const vtkStdString& name, int flag);
461 
463  {
464  this->FastPathObjectType = type;
465  }
466  void SetFastPathObjectId(vtkIdType id) { this->FastPathObjectId = id; }
467  vtkSetStringMacro(FastPathIdType);
468 
470 
479 
488 
495  void SetInitialObjectStatus(int otype, const char* name, int stat);
496 
502  void SetInitialObjectArrayStatus(int otype, const char* name, int stat);
503 
505 
507 
508 protected:
511 
513  void BuildSIL();
514 
518  int nn, char** np, vtksys::RegularExpression& re, vtkStdString& field, vtkStdString& ele);
519 
521  void GlomArrayNames(int i, int num_obj, int num_vars, char** var_names, int* truth_tab);
522 
525 
542  int AssembleOutputConnectivity(vtkIdType timeStep, int otyp, int oidx, int conntypidx,
543  BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
551  vtkIdType timeStep, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
556  vtkIdType timeStep, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
561  vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
567  vtkIdType timeStep, int otyp, int oidx, vtkUnstructuredGrid* output);
570  vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
578  vtkIdType timeStep, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
580  vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
584 
602 
605 
608  BlockInfoType* binfo, vtkIntArray* facesPerCell, vtkIdTypeArray* exoCellConn);
609 
611  void InsertBlockCells(int otyp, int obj, int conn_type, int timeStep, BlockInfoType* binfop);
612 
614  void InsertSetCells(int otyp, int obj, int conn_type, int timeStep, SetInfoType* sinfop);
615 
618 
620  void InsertSetNodeCopies(vtkIdTypeArray* refs, int otyp, int obj, SetInfoType* sinfo);
621 
623  void InsertSetCellCopies(vtkIdTypeArray* refs, int otyp, int obj, SetInfoType* sinfo);
624 
626  void InsertSetSides(vtkIdTypeArray* refs, int otyp, int obj, SetInfoType* sinfo);
627 
634 
640 
646 
652  int GetNumberOfObjectsAtTypeIndex(int typeIndex);
653 
661  ObjectInfoType* GetObjectInfo(int typeIndex, int objectIndex);
662 
669  ObjectInfoType* GetSortedObjectInfo(int objectType, int objectIndex);
670 
677  ObjectInfoType* GetUnsortedObjectInfo(int objectType, int objectIndex);
678 
683  int GetBlockIndexFromFileGlobalId(int otyp, int refId);
684 
690 
695 
698 
702  ArrayInfoType* FindArrayInfoByName(int otyp, const char* name);
703 
707  int IsObjectTypeBlock(int otyp);
708  int IsObjectTypeSet(int otyp);
709  int IsObjectTypeMap(int otyp);
710 
717 
722 
727 
733  void RemoveBeginningAndTrailingSpaces(int len, char** names, int maxNameLength);
734 
737 
741  std::map<int, std::vector<BlockInfoType>> BlockInfo;
745  std::map<int, std::vector<SetInfoType>> SetInfo;
751  std::map<int, std::vector<MapInfoType>> MapInfo;
752 
753  std::vector<PartInfoType> PartInfo;
754  std::vector<MaterialInfoType> MaterialInfo;
755  std::vector<AssemblyInfoType> AssemblyInfo;
756 
761  std::map<int, std::vector<int>> SortedObjectIndices;
763  // defined on that type.
764  std::map<int, std::vector<ArrayInfoType>> ArrayInfo;
765 
770  std::map<int, std::vector<ArrayInfoType>> InitialArrayInfo;
771 
776  std::map<int, std::vector<ObjectInfoType>> InitialObjectInfo;
777 
781 
786 
788  int Exoid;
789 
791  struct ex_init_params ModelParameters;
792 
794  std::vector<double> Times;
796 
801 
809 
813  int FileId;
814 
817  //
819  double CacheSize;
820 
825 
827 
840 
845 
847 
856  std::map<int, std::vector<std::vector<vtkIdType>>> PolyhedralFaceConnArrays;
857 
861 
863 
864 private:
866  void operator=(const vtkExodusIIReaderPrivate&) = delete;
867 };
868 
869 #endif
870 #endif
871 #endif // vtkExodusIIReaderPrivate_h
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
internal parser used by vtkExodusIIReader.
This class holds metadata for an Exodus file.
int AssembleOutputPointMaps(vtkIdType timeStep, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Add maps to an output mesh.
const char * GetObjectAttributeName(int objectType, int objectIndex, int attributeIndex)
void SetFastPathObjectId(vtkIdType id)
void SetAssemblyStatus(const vtkStdString &name, int flag)
vtkTimeStamp InformationTimeStamp
Time stamp from last time we were in RequestInformation.
vtkExodusIIReader::ObjectType FastPathObjectType
static const char * GetGlobalVariableValuesArrayName()
void Reset()
Reset the class so that another file may be read.
virtual void SetDisplacementMagnitude(double s)
void InsertSetCellCopies(vtkIdTypeArray *refs, int otyp, int obj, SetInfoType *sinfo)
Insert cells referenced by an edge, face, or element set.
void ResetCache()
Clears out any data in the cache and restores it to its initial state.
const char * GetAssemblyName(int idx)
static const char * GetFileIdArrayName()
int IsObjectTypeSet(int otyp)
int CloseFile()
Close any ExodusII file currently open for reading. Returns 0 on success.
ObjectInfoType * GetSortedObjectInfo(int objectType, int objectIndex)
Return a pointer to the ObjectInfo of the specified type and index, but using indices sorted by objec...
void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status)
vtkDataArray * FindDisplacementVectors(int timeStep)
void AddPointArray(vtkDataArray *src, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Add a point array to an output grid's point data, squeezing if necessary.
ArraySourceTypes
Tags to indicate the source of values for an array.
std::vector< double > Times
A list of time steps for which results variables are stored.
int GetObjectArrayStatus(int otype, int i)
For a given object type, returns the status of the i-th array.
int OpenFile(const char *filename)
Open an ExodusII file for reading. Returns 0 on success.
int GetAssemblyStatus(int idx)
int GetMaterialStatus(const vtkStdString &name)
int AssembleOutputPoints(vtkIdType timeStep, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Fill the output grid's point coordinates array.
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.
void PrepareGeneratedArrayInfo()
Add generated array information to array info lists.
void SetCacheSize(double size)
Set the size of the cache in MiB.
const char * GetObjectName(int otype, int i)
For a given object type, returns the name of the i-th object.
void SetMaterialStatus(const vtkStdString &name, int flag)
vtkDataArray * GetCacheOrRead(vtkExodusIICacheKey)
Return an array for the specified cache key.
vtkExodusIIReader * Parent
Pointer to owning reader...
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.
int RequestData(vtkIdType timeStep, vtkMultiBlockDataSet *output)
Read requested data and store in unstructured grid.
int GetObjectTypeIndexFromObjectType(int otyp)
Return the index of an object type (in a private list of all object types).
int GetBlockConnTypeFromBlockType(int btyp)
Given a block type (EDGE_BLOCK, ...), return the associated block connectivity type (EDGE_BLOCK_CONN,...
BlockInfoType * GetBlockFromFileGlobalId(int otyp, int refId)
Get the block containing the entity referenced by the specified file-global ID.
static const char * GetObjectIdArrayName()
int GetUnsortedObjectStatus(int otype, int i)
For a given object type, returns the status of the i-th object, where i is an index into the unsorted...
int GetObjectTypeFromMapType(int mtyp)
Given a map type (NODE_MAP, EDGE_MAP, ...) return the associated object type (NODAL,...
static const char * GetImplicitNodeIdArrayName()
void SetPartStatus(int idx, int on)
static const char * GetGlobalNodeIdArrayName()
vtkIdType GetPolyhedronFaceConnectivity(vtkIdType fileLocalFaceId, vtkIdType *&facePtIds)
Fetch the face-connectivity for one face of one polyhedron.
int GetBlockIndexFromFileGlobalId(int otyp, int refId)
Get the index of the block containing the entity referenced by the specified file-global ID.
static const char * GetGlobalVariableNamesArrayName()
int IsObjectTypeBlock(int otyp)
Does the specified object type match? Avoid using these...
void FreePolyhedronFaceArrays()
Free any arrays held by PolyhedralFaceConnArrays (for polyhedral-face-connectivity lookup).
int GetNumberOfObjectAttributes(int objectType, int objectIndex)
Unlike object arrays, attributes are only defined over blocks (not sets) and are defined on a per-blo...
int GetNumberOfObjectArrayComponents(int otype, int i)
For a given object type, returns the number of components of the i-th array.
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 GetMapTypeFromObjectType(int otyp)
void SetAssemblyStatus(int idx, int on)
int GetObjectAttributeIndex(int objectType, int objectIndex, const char *attribName)
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 InsertBlockCells(int otyp, int obj, int conn_type, int timeStep, BlockInfoType *binfop)
Insert cells from a specified block into a mesh.
int GetObjectStatus(int otype, int i)
For a given object type, returns the status of the i-th object.
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.
double ModeShapeTime
The time value.
ArrayInfoType * FindArrayInfoByName(int otyp, const char *name)
Find an ArrayInfo object for a specific object type using the name as a key.
void DetermineVtkCellType(BlockInfoType &binfo)
Determine the VTK cell type for a given edge/face/element block.
ObjectInfoType * GetObjectInfo(int typeIndex, int objectIndex)
Return a pointer to the ObjectInfo of the specified type and index.
int AssembleArraysOverTime(vtkMultiBlockDataSet *output)
Add fast-path time-varying data to field data of an output block or set.
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::vector< AssemblyInfoType > AssemblyInfo
int GetObjectId(int otype, int i)
For a given object type, return the user-assigned ID of the i-th object.
void BuildSIL()
Build SIL. This must be called only after RequestInformation().
void InsertSetCells(int otyp, int obj, int conn_type, int timeStep, SetInfoType *sinfop)
Insert cells from a specified set into a mesh.
int IsObjectTypeMap(int otyp)
void SetSqueezePoints(int sp)
Set whether subsequent RequestData() calls will produce the minimal point set required to represent t...
int AssembleOutputPointArrays(vtkIdType timeStep, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Add the requested arrays to the output grid's point data.
int GetNumberOfNodes()
Return the number of nodes in the output (depends on SqueezePoints)
const char * GetPartName(int idx)
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 ...
float ExodusVersion
The version of Exodus that wrote the currently open file (or a negative number otherwise).
int GetConnTypeIndexFromConnType(int ctyp)
Return the index of an object type (in a private list of all object types).
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 SetFastPathObjectType(vtkExodusIIReader::ObjectType type)
const char * GetMaterialName(int idx)
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)
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.
int GetAssemblyStatus(const vtkStdString &name)
std::map< int, std::vector< std::vector< vtkIdType > > > PolyhedralFaceConnArrays
Face connectivity for polyhedra.
vtkMutableDirectedGraph * SIL
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 RemoveBeginningAndTrailingSpaces(int len, char **names, int maxNameLength)
Function to trim space from names retrieved with ex_get_var_names.
std::vector< MaterialInfoType > MaterialInfo
int GetNumberOfObjectsOfType(int otype)
Returns the number of objects of a given type (e.g., EX_ELEM_BLOCK, EX_NODE_SET, ....
static const char * GetGlobalElementIdArrayName()
int GetObjectSize(int otype, int i)
For a given object type, return the size of the i-th object.
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.
std::map< int, std::vector< ArrayInfoType > > ArrayInfo
Maps an object type (EX_ELEM_BLOCK, EX_NODE_SET, ...) to a list of arrays.
void SetPartStatus(const vtkStdString &name, int flag)
void SetObjectArrayStatus(int otype, int i, int stat)
For a given object type, sets the status of the i-th array.
static const char * GetImplicitElementIdArrayName()
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetPartStatus(const vtkStdString &name)
~vtkExodusIIReaderPrivate() override
vtkMutableDirectedGraph * GetSIL()
Returns the SIL. This valid only after BuildSIL() has been called.
int RequestInformation()
Get metadata for an open file with handle exoid.
void ClearConnectivityCaches()
Delete any cached connectivity information (for all blocks and sets)
int AppWordSize
These aren't the variables you're looking for.
void ResetSettings()
Return user-specified variables to their default values.
const char * GetObjectArrayName(int otype, int i)
For a given object type, returns the name of the i-th array.
void InsertSetSides(vtkIdTypeArray *refs, int otyp, int obj, SetInfoType *sinfo)
Insert cells referenced by a side set.
int GetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex)
void SetUnsortedObjectStatus(int otype, int i, int stat)
For a given object type, sets the status of the i-th object, where i is an index into the unsorted ob...
int GetNumberOfObjectArraysOfType(int otype)
Returns the number of arrays defined over objects of a given type (e.g., EX_ELEM_BLOCK,...
int SqueezePoints
Should the reader output only points used by elements in the output mesh, or all the points.
int AssembleOutputCellMaps(vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
int AssembleOutputCellArrays(vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType *bsinfop, vtkUnstructuredGrid *output)
Add the requested arrays to the output grid's cell data.
const struct ex_init_params * GetModelParams() const
void InsertBlockPolyhedra(BlockInfoType *binfo, vtkIntArray *facesPerCell, vtkIdTypeArray *exoCellConn)
Insert polyhedral cells (called from InsertBlockCells when a block is polyhedral).
vtkExodusIIReaderParser * Parser
int GetPartStatus(int idx)
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.
double CacheSize
The size of the cache in MiB.
int AssembleOutputProceduralArrays(vtkIdType timeStep, int otyp, int oidx, vtkUnstructuredGrid *output)
Add procedurally generated arrays to an output mesh.
const char * GetPartBlockInfo(int idx)
void SetMaterialStatus(int idx, int on)
virtual void SetParser(vtkExodusIIReaderParser *)
void InsertSetNodeCopies(vtkIdTypeArray *refs, int otyp, int obj, SetInfoType *sinfo)
Insert cells referenced by a node set.
int GetNumberOfTimeSteps()
Return the number of time steps in the open file.
int SetUpEmptyGrid(vtkMultiBlockDataSet *output)
Description: Prepare a data set with the proper structure and arrays but no cells.
int GetMaterialStatus(int idx)
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.
vtkIdType GetSqueezePointId(BlockSetInfoType *bsinfop, int i)
Find or create a new SqueezePoint ID (unique sequential list of points referenced by cells in blocks/...
std::vector< PartInfoType > PartInfo
virtual void SetApplyDisplacements(vtkTypeBool d)
vtkExodusIICache * Cache
A least-recently-used cache to hold raw arrays.
void SetTimesOverrides(const std::vector< double > &times)
GlomTypes
Tags to indicate how single-component Exodus arrays are glommed (aggregated) into multi-component VTK...
int GetSetTypeFromSetConnType(int sctyp)
Given a set connectivity type (NODE_SET_CONN, ...), return the associated object type (NODE_SET,...
int Exoid
The handle of the currently open file.
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 GetNumberOfObjectsAtTypeIndex(int typeIndex)
Return the number of objects of the given type.
void SetObjectStatus(int otype, int i, int stat)
For a given object type, sets the status of the i-th object.
static vtkExodusIIReaderPrivate * New()
Read exodus 2 files .ex2.
ObjectType
Extra cell data array that can be generated.
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:34
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:40
Composite dataset that organizes datasets into blocks.
An editable directed graph.
abstract base class for most VTK objects
Definition: vtkObject.h:63
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
Read Exodus II files (.exii)
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:36
record modification and/or execution time
Definition: vtkTimeStamp.h:33
dataset represents arbitrary combinations of all possible cell types
@ field
Definition: vtkX3D.h:183
@ on
Definition: vtkX3D.h:445
@ info
Definition: vtkX3D.h:382
@ type
Definition: vtkX3D.h:522
@ name
Definition: vtkX3D.h:225
@ size
Definition: vtkX3D.h:259
A struct to hold information about time-varying arrays.
int StorageType
Storage type of array (a type that can be passed to vtkDataArray::Create())
void Reset()
Clear all the structure members.
int GlomType
The type of "glomming" performed.
std::vector< int > ObjectTruth
A map describing which objects the variable is defined on.
vtkStdString Name
The name of the array.
int Source
The source of the array (Result or Attribute)
std::vector< int > OriginalIndices
The index of each component of the array as ordered by the Exodus file.
int Status
Whether or not the array should be loaded by RequestData.
std::vector< vtkStdString > OriginalNames
The name of each component of the array as defined by the Exodus file.
int Components
The number of components in the array.
A struct to hold information about Exodus blocks.
A struct to hold information about Exodus blocks or sets (they have some members in common)
std::map< vtkIdType, vtkIdType > PointMap
A map from nodal IDs in an Exodus file to nodal IDs in the output mesh.
vtkIdType FileOffset
Id (1-based) of first entry in file-local list across all blocks in file.
BlockSetInfoType & operator=(const BlockSetInfoType &block)
std::map< vtkIdType, vtkIdType > ReversePointMap
A map from nodal ids in the output mesh to those in an Exodus file.
vtkUnstructuredGrid * CachedConnectivity
Cached cell connectivity arrays for mesh.
BlockSetInfoType(const BlockSetInfoType &block)
vtkIdType NextSqueezePoint
The next vtk ID to use for a connectivity entry when point squeezing is on and no point ID exists.
A struct to hold information about Exodus maps.
A struct to hold information about Exodus objects (blocks, sets, maps)
int Size
Number of entries in this block.
int Id
User-assigned identification number.
int Status
Should the reader load this block?
A struct to hold information about Exodus blocks.
A struct to hold information about Exodus sets.
int vtkTypeBool
Definition: vtkABI.h:69
#define VTK_DEPRECATED_IN_9_1_0(reason)
int vtkIdType
Definition: vtkType.h:332