VTK  9.4.20250110
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
20VTK_ABI_NAMESPACE_BEGIN
21class vtkDataArray;
23class vtkIdTypeArray;
26class vtkTypeInt64Array;
28
32class VTKIOEXODUS_EXPORT vtkExodusIIReaderPrivate : public vtkObject
33{
34public:
36 void PrintSelf(ostream& os, vtkIndent indent) override;
38 // virtual void Modified();
39
41 int OpenFile(const char* filename);
42
44 int CloseFile();
45
48
50 vtkMutableDirectedGraph* GetSIL() { return this->SIL; }
51
54
60
72 void Reset();
73
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
109
115
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
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
265
266 const struct ex_init_params* GetModelParams() const { return &this->ModelParameters; }
267
269 struct VTKIOEXODUS_EXPORT ArrayInfoType
270 {
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;
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; }
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...)
368 // Number of points per cell as used by VTK
369 // -- not what's in the file (i.e., BdsPerEntry[0] >= PointsPerCell)
371 };
372
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
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
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
448 const char* GetMaterialName(int idx);
449 int GetMaterialStatus(int idx);
451 void SetMaterialStatus(int idx, int on);
452 void SetMaterialStatus(const vtkStdString& name, int flag);
453
455 const char* GetAssemblyName(int idx);
456 int GetAssemblyStatus(int idx);
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
469
478
487
494 void SetInitialObjectStatus(int otype, const char* name, int stat);
495
501 void SetInitialObjectArrayStatus(int otype, const char* name, int stat);
502
504
506
507protected:
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
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
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
639
645
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
689
694
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
716
721
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
813
816 //
818 double CacheSize;
819
824
826
839
844
846
855 std::map<int, std::vector<std::vector<vtkIdType>>> PolyhedralFaceConnArrays;
856
860
862
863private:
865 void operator=(const vtkExodusIIReaderPrivate&) = delete;
866};
867
868VTK_ABI_NAMESPACE_END
869#endif // vtkExodusIIReaderPrivate_h
abstract superclass for arrays of numeric data
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.
static const char * GetImplicitNodeIdArrayName()
void SetAssemblyStatus(const vtkStdString &name, int flag)
vtkTimeStamp InformationTimeStamp
Time stamp from last time we were in RequestInformation.
vtkExodusIIReader::ObjectType FastPathObjectType
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.
int IsObjectTypeSet(int otyp)
int CloseFile()
Close any ExodusII file currently open for reading. Returns 0 on success.
void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status)
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.
static const char * GetGlobalNodeIdArrayName()
void PrepareGeneratedArrayInfo()
Add generated array information to array info lists.
void SetCacheSize(double size)
Set the size of the cache in MiB.
void SetMaterialStatus(const vtkStdString &name, int flag)
vtkExodusIIReader * Parent
Pointer to owning reader... this is not registered in order to avoid circular references.
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,...
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,...
void SetPartStatus(int idx, int on)
const char * GetAssemblyName(int idx)
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.
int IsObjectTypeBlock(int otyp)
Does the specified object type match? Avoid using these... they aren't robust against new types being...
void FreePolyhedronFaceArrays()
Free any arrays held by PolyhedralFaceConnArrays (for polyhedral-face-connectivity lookup).
const struct ex_init_params * GetModelParams() const
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.
const char * GetPartBlockInfo(int idx)
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.
static vtkExodusIIReaderPrivate * New()
void DetermineVtkCellType(BlockInfoType &binfo)
Determine the VTK cell type for a given edge/face/element block.
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...
static const char * GetObjectIdArrayName()
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)
vtkMutableDirectedGraph * GetSIL()
Returns the SIL. This valid only after BuildSIL() has been called.
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)
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)
BlockInfoType * GetBlockFromFileGlobalId(int otyp, int refId)
Get the block containing the entity referenced by the specified file-global ID.
std::map< int, std::vector< std::vector< vtkIdType > > > PolyhedralFaceConnArrays
Face connectivity for polyhedra.
ObjectInfoType * GetUnsortedObjectInfo(int objectType, int objectIndex)
Return a pointer to the ObjectInfo of the specified type and index, but using indices sorted by objec...
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, ....
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.
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
int RequestInformation()
Get metadata for an open file with handle exoid.
void ClearConnectivityCaches()
Delete any cached connectivity information (for all blocks and sets)
static const char * GetFileIdArrayName()
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...
const char * GetObjectName(int otype, int i)
For a given object type, returns the name of the i-th object.
int GetNumberOfObjectArraysOfType(int otype)
Returns the number of arrays defined over objects of a given type (e.g., EX_ELEM_BLOCK,...
ObjectInfoType * GetObjectInfo(int typeIndex, int objectIndex)
Return a pointer to the ObjectInfo of the specified type and index.
int SqueezePoints
Should the reader output only points used by elements in the output mesh, or all the points.
static const char * GetGlobalVariableNamesArrayName()
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.
void InsertBlockPolyhedra(BlockInfoType *binfo, vtkIntArray *facesPerCell, vtkIdTypeArray *exoCellConn)
Insert polyhedral cells (called from InsertBlockCells when a block is polyhedral).
static const char * GetGlobalElementIdArrayName()
ObjectInfoType * GetSortedObjectInfo(int objectType, int objectIndex)
Return a pointer to the ObjectInfo of the specified type and index, but using indices sorted by objec...
vtkExodusIIReaderParser * Parser
static const char * GetImplicitElementIdArrayName()
int GetPartStatus(int idx)
const char * GetPartName(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.
void SetMaterialStatus(int idx, int on)
static const char * GetGlobalVariableValuesArrayName()
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.
vtkDataArray * FindDisplacementVectors(int timeStep)
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
ArrayInfoType * FindArrayInfoByName(int otyp, const char *name)
Find an ArrayInfo object for a specific object type using the name as a key.
virtual void SetApplyDisplacements(vtkTypeBool d)
vtkDataArray * GetCacheOrRead(vtkExodusIICacheKey)
Return an array for the specified cache key.
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.
const char * GetObjectAttributeName(int objectType, int objectIndex, int attributeIndex)
void SetObjectStatus(int otype, int i, int stat)
For a given object type, sets the status of the i-th object.
const char * GetMaterialName(int idx)
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:108
dynamic, self-adjusting array of int
Composite dataset that organizes datasets into blocks.
An editable directed graph.
abstract base class for most VTK objects
Definition vtkObject.h:162
Read Exodus II files (.exii)
Wrapper around std::string to keep symbols short.
record modification and/or execution time
dataset represents arbitrary combinations of all possible cell types
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.
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.
BlockSetInfoType & operator=(const BlockSetInfoType &block)
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:64
int vtkIdType
Definition vtkType.h:315