VTK  9.3.1
vtkExodusIIReaderPrivate.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-License-Identifier: BSD-3-Clause
3 #ifndef vtkExodusIIReaderPrivate_h
4 #define vtkExodusIIReaderPrivate_h
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" // for vtkExodusIICacheKey
10 #include "vtkExodusIIReader.h" // for vtkExodusIIReader
11 #include "vtkObject.h"
12 #include "vtkStdString.h" // for vtkStdString
13 #include "vtksys/RegularExpression.hxx" // for vtksys::RegularExpression
14 
15 #include <map> // for std::map
16 #include <vector> // for std::vector
17 
18 #include "vtkIOExodusModule.h" // For export macro
19 #include "vtk_exodusII.h" // for exodus APIs
20 VTK_ABI_NAMESPACE_BEGIN
21 class vtkDataArray;
23 class vtkIdTypeArray;
26 class vtkTypeInt64Array;
28 
32 class VTKIOEXODUS_EXPORT vtkExodusIIReaderPrivate : public vtkObject
33 {
34 public:
35  static vtkExodusIIReaderPrivate* New();
36  void PrintSelf(ostream& os, vtkIndent indent) override;
38  // virtual void Modified();
39 
41  int OpenFile(const char* filename);
42 
44  int CloseFile();
45 
47  int RequestInformation();
48 
50  vtkMutableDirectedGraph* GetSIL() { return this->SIL; }
51 
53  int RequestData(vtkIdType timeStep, vtkMultiBlockDataSet* output);
54 
59  int SetUpEmptyGrid(vtkMultiBlockDataSet* output);
60 
72  void Reset();
73 
78  void ResetSettings();
79 
81  void ResetCache();
82 
84  void SetCacheSize(double size);
85 
87  vtkGetMacro(CacheSize, double);
88 
93  int GetNumberOfTimeSteps() { return (int)this->Times.size(); }
94 
97  vtkGetMacro(SqueezePoints, int);
98 
101  void SetSqueezePoints(int sp);
102 
105  vtkBooleanMacro(SqueezePoints, int);
106 
108  int GetNumberOfNodes();
109 
114  int GetNumberOfObjectsOfType(int otype);
115 
126  int GetNumberOfObjectArraysOfType(int otype);
127 
132  const char* GetObjectName(int otype, int i);
133  using Superclass::GetObjectName;
134 
139  int GetObjectId(int otype, int i);
140 
147  int GetObjectSize(int otype, int i);
148 
153  int GetObjectStatus(int otype, int i);
154 
160  int GetUnsortedObjectStatus(int otype, int i);
161 
166  void SetObjectStatus(int otype, int i, int stat);
167 
173  void SetUnsortedObjectStatus(int otype, int i, int stat);
174 
179  const char* GetObjectArrayName(int otype, int i);
180 
185  int GetNumberOfObjectArrayComponents(int otype, int i);
186 
191  int GetObjectArrayStatus(int otype, int i);
192 
197  void SetObjectArrayStatus(int otype, int i, int stat);
198 
205  int GetNumberOfObjectAttributes(int objectType, int objectIndex);
206  const char* GetObjectAttributeName(int objectType, int objectIndex, int attributeIndex);
207  int GetObjectAttributeIndex(int objectType, int objectIndex, const char* attribName);
208  int GetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex);
209  void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status);
210 
212  vtkGetMacro(GenerateObjectIdArray, vtkTypeBool);
213  vtkSetMacro(GenerateObjectIdArray, vtkTypeBool);
214  static const char* GetObjectIdArrayName() { return "ObjectId"; }
215 
216  vtkSetMacro(GenerateGlobalElementIdArray, vtkTypeBool);
217  vtkGetMacro(GenerateGlobalElementIdArray, vtkTypeBool);
218  static const char* GetGlobalElementIdArrayName() { return "GlobalElementId"; }
219 
220  vtkSetMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
221  vtkGetMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
222  static const char* GetGlobalNodeIdArrayName() { return "GlobalNodeId"; }
223 
224  vtkSetMacro(GenerateImplicitElementIdArray, vtkTypeBool);
225  vtkGetMacro(GenerateImplicitElementIdArray, vtkTypeBool);
226  static const char* GetImplicitElementIdArrayName() { return "ImplicitElementId"; }
227 
228  vtkSetMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
229  vtkGetMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
230  static const char* GetImplicitNodeIdArrayName() { return "ImplicitNodeId"; }
231 
235  vtkSetMacro(GenerateFileIdArray, vtkTypeBool);
236  vtkGetMacro(GenerateFileIdArray, vtkTypeBool);
237  static const char* GetFileIdArrayName() { return "FileId"; }
238 
240  vtkSetMacro(FileId, int);
241  vtkGetMacro(FileId, int);
242 
243  static const char* GetGlobalVariableValuesArrayName() { return "GlobalVariableValues"; }
244  static const char* GetGlobalVariableNamesArrayName() { return "GlobalVariableNames"; }
245 
246  virtual void SetApplyDisplacements(vtkTypeBool d);
247  vtkGetMacro(ApplyDisplacements, vtkTypeBool);
248 
249  virtual void SetDisplacementMagnitude(double s);
250  vtkGetMacro(DisplacementMagnitude, double);
251 
252  vtkSetMacro(HasModeShapes, int);
253  vtkGetMacro(HasModeShapes, int);
254 
255  vtkSetMacro(ModeShapeTime, double);
256  vtkGetMacro(ModeShapeTime, double);
257 
258  vtkSetMacro(AnimateModeShapes, int);
259  vtkGetMacro(AnimateModeShapes, int);
260 
261  vtkSetMacro(IgnoreFileTime, bool);
262  vtkGetMacro(IgnoreFileTime, bool);
263 
264  vtkDataArray* FindDisplacementVectors(int timeStep);
265 
266  const struct ex_init_params* GetModelParams() const { return &this->ModelParameters; }
267 
269  struct VTKIOEXODUS_EXPORT ArrayInfoType
270  {
281  int GlomType;
286  int Source;
288  int Status;
291  std::vector<vtkStdString> OriginalNames;
294  std::vector<int> OriginalIndices;
303  std::vector<int> ObjectTruth;
305  void Reset();
306  };
307 
309  struct VTKIOEXODUS_EXPORT ObjectInfoType
310  {
312  int Size;
314  int Status;
316  int Id;
319  };
320 
322  struct VTKIOEXODUS_EXPORT MapInfoType : public ObjectInfoType
323  {
324  };
325 
328  struct VTKIOEXODUS_EXPORT BlockSetInfoType : public ObjectInfoType
329  {
336  std::map<vtkIdType, vtkIdType> PointMap;
341  std::map<vtkIdType, vtkIdType> ReversePointMap;
348 
349  BlockSetInfoType() { this->CachedConnectivity = nullptr; }
350  BlockSetInfoType(const BlockSetInfoType& block);
351  ~BlockSetInfoType();
352  BlockSetInfoType& operator=(const BlockSetInfoType& block);
353  };
354 
356  struct VTKIOEXODUS_EXPORT BlockInfoType : public BlockSetInfoType
357  {
358  vtkStdString OriginalName; // useful to reset the name if XML metadata is invalid.
360  // number of boundaries per entry
361  // The index is the dimensionality of the entry. 0=node, 1=edge, 2=face
362  int64_t BdsPerEntry[3];
364  std::vector<vtkStdString> AttributeNames;
365  std::vector<int> AttributeStatus;
366  // VTK cell type (a function of TypeName and BdsPerEntry...)
367  int CellType;
368  // Number of points per cell as used by VTK
369  // -- not what's in the file (i.e., BdsPerEntry[0] >= PointsPerCell)
371  };
372 
374  struct PartInfoType : public ObjectInfoType
375  {
376  std::vector<int> BlockIndices;
377  };
379  {
380  std::vector<int> BlockIndices;
381  };
383  {
384  std::vector<int> BlockIndices;
385  };
386 
389  {
390  int DistFact; // Number of distribution factors
391  // (for the entire block, not per array or entry)
392  };
393 
397  {
398  Scalar = 0,
399  Vector2 = 1,
400  Vector3 = 2,
401  SymmetricTensor = 3,
402  // (order xx, yy, zz, xy, yz, zx)
403  IntegrationPoint = 4
404  };
405 
408  {
409  Result = 0,
410  // (that vary over time)
411  Attribute = 1,
412  // (constants over time)
413  Map = 2,
414  Generated = 3
415  };
416 
419 
420  friend class vtkExodusIIReader;
421  friend class vtkPExodusIIReader;
422 
423  virtual void SetParser(vtkExodusIIReaderParser*);
424  vtkGetObjectMacro(Parser, vtkExodusIIReaderParser);
425 
426  // BUG #15632: This method allows vtkPExodusIIReader to pass time information
427  // from one spatial file to another and avoiding have to read it for each of
428  // the files.
429  void SetTimesOverrides(const std::vector<double>& times)
430  {
431  this->Times = times;
432  this->SkipUpdateTimeInformation = true;
433  }
434 
435  // Because Parts, Materials, and assemblies are not stored as arrays,
436  // but rather as maps to the element blocks they make up,
437  // we cannot use the Get|SetObject__() methods directly.
438 
439  int GetNumberOfParts();
440  const char* GetPartName(int idx);
441  const char* GetPartBlockInfo(int idx);
442  int GetPartStatus(int idx);
443  int GetPartStatus(const vtkStdString& name);
444  void SetPartStatus(int idx, int on);
445  void SetPartStatus(const vtkStdString& name, int flag);
446 
447  int GetNumberOfMaterials();
448  const char* GetMaterialName(int idx);
449  int GetMaterialStatus(int idx);
450  int GetMaterialStatus(const vtkStdString& name);
451  void SetMaterialStatus(int idx, int on);
452  void SetMaterialStatus(const vtkStdString& name, int flag);
453 
454  int GetNumberOfAssemblies();
455  const char* GetAssemblyName(int idx);
456  int GetAssemblyStatus(int idx);
457  int GetAssemblyStatus(const vtkStdString& name);
458  void SetAssemblyStatus(int idx, int on);
459  void SetAssemblyStatus(const vtkStdString& name, int flag);
460 
462  {
463  this->FastPathObjectType = type;
464  }
465  void SetFastPathObjectId(vtkIdType id) { this->FastPathObjectId = id; }
466  vtkSetStringMacro(FastPathIdType);
467 
468  bool IsXMLMetadataValid();
469 
477  void GetInitialObjectStatus(int otype, ObjectInfoType* info);
478 
486  void GetInitialObjectArrayStatus(int otype, ArrayInfoType* info);
487 
494  void SetInitialObjectStatus(int otype, const char* name, int stat);
495 
501  void SetInitialObjectArrayStatus(int otype, const char* name, int stat);
502 
503  int UpdateTimeInformation();
504 
506 
507 protected:
509  ~vtkExodusIIReaderPrivate() override;
510 
512  void BuildSIL();
513 
517  int nn, char** np, vtksys::RegularExpression& re, vtkStdString& field, vtkStdString& ele);
518 
520  void GlomArrayNames(int i, int num_obj, int num_vars, char** var_names, int* truth_tab);
521 
524 
541  int AssembleOutputConnectivity(vtkIdType timeStep, int otyp, int oidx, int conntypidx,
542  BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
550  vtkIdType timeStep, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
555  vtkIdType timeStep, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
560  vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
566  vtkIdType timeStep, int otyp, int oidx, vtkUnstructuredGrid* output);
569  vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
577  vtkIdType timeStep, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
579  vtkIdType timeStep, int otyp, int oidx, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
583 
600  vtkIdType GetPolyhedronFaceConnectivity(vtkIdType fileLocalFaceId, vtkIdType*& facePtIds);
601 
604 
607  BlockInfoType* binfo, vtkIntArray* facesPerCell, vtkIdTypeArray* exoCellConn);
608 
610  void InsertBlockCells(int otyp, int obj, int conn_type, int timeStep, BlockInfoType* binfop);
611 
613  void InsertSetCells(int otyp, int obj, int conn_type, int timeStep, SetInfoType* sinfop);
614 
616  void AddPointArray(vtkDataArray* src, BlockSetInfoType* bsinfop, vtkUnstructuredGrid* output);
617 
619  void InsertSetNodeCopies(vtkIdTypeArray* refs, int otyp, int obj, SetInfoType* sinfo);
620 
622  void InsertSetCellCopies(vtkIdTypeArray* refs, int otyp, int obj, SetInfoType* sinfo);
623 
625  void InsertSetSides(vtkIdTypeArray* refs, int otyp, int obj, SetInfoType* sinfo);
626 
633 
638  int GetConnTypeIndexFromConnType(int ctyp);
639 
644  int GetObjectTypeIndexFromObjectType(int otyp);
645 
651  int GetNumberOfObjectsAtTypeIndex(int typeIndex);
652 
660  ObjectInfoType* GetObjectInfo(int typeIndex, int objectIndex);
661 
668  ObjectInfoType* GetSortedObjectInfo(int objectType, int objectIndex);
669 
676  ObjectInfoType* GetUnsortedObjectInfo(int objectType, int objectIndex);
677 
682  int GetBlockIndexFromFileGlobalId(int otyp, int refId);
683 
688  BlockInfoType* GetBlockFromFileGlobalId(int otyp, int refId);
689 
694 
696  void DetermineVtkCellType(BlockInfoType& binfo);
697 
701  ArrayInfoType* FindArrayInfoByName(int otyp, const char* name);
702 
706  int IsObjectTypeBlock(int otyp);
707  int IsObjectTypeSet(int otyp);
708  int IsObjectTypeMap(int otyp);
709 
713  int GetObjectTypeFromMapType(int mtyp);
714  int GetMapTypeFromObjectType(int otyp);
715  int GetTemporalTypeFromObjectType(int otyp);
716 
720  int GetSetTypeFromSetConnType(int sctyp);
721 
725  int GetBlockConnTypeFromBlockType(int btyp);
726 
732  void RemoveBeginningAndTrailingSpaces(int len, char** names, int maxNameLength);
733 
736 
740  std::map<int, std::vector<BlockInfoType>> BlockInfo;
744  std::map<int, std::vector<SetInfoType>> SetInfo;
750  std::map<int, std::vector<MapInfoType>> MapInfo;
751 
752  std::vector<PartInfoType> PartInfo;
753  std::vector<MaterialInfoType> MaterialInfo;
754  std::vector<AssemblyInfoType> AssemblyInfo;
755 
760  std::map<int, std::vector<int>> SortedObjectIndices;
762  // defined on that type.
763  std::map<int, std::vector<ArrayInfoType>> ArrayInfo;
764 
769  std::map<int, std::vector<ArrayInfoType>> InitialArrayInfo;
770 
775  std::map<int, std::vector<ObjectInfoType>> InitialObjectInfo;
776 
780 
785 
787  int Exoid;
788 
790  struct ex_init_params ModelParameters;
791 
793  std::vector<double> Times;
795 
800 
808 
812  int FileId;
813 
816  //
818  double CacheSize;
819 
824 
826 
839 
844 
846 
855  std::map<int, std::vector<std::vector<vtkIdType>>> PolyhedralFaceConnArrays;
856 
860 
862 
863 private:
865  void operator=(const vtkExodusIIReaderPrivate&) = delete;
866 };
867 
868 VTK_ABI_NAMESPACE_END
869 #endif // 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:28
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:51
A struct to hold information about Exodus objects (blocks, sets, maps)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
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'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:24
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'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:315
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:64
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:34
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:28
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:44
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.
const struct ex_init_params * GetModelParams() const
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) ...
virtual std::string GetObjectName() const
Set/get the name of this object for reporting purposes.
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'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)
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't the variables you'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'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)