VTK  9.3.0
vtkExodusIIReader.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3// SPDX-License-Identifier: BSD-3-Clause
4
34#ifndef vtkExodusIIReader_h
35#define vtkExodusIIReader_h
36
37#include "vtkIOExodusModule.h" // For export macro
39
40VTK_ABI_NAMESPACE_BEGIN
41class vtkDataArray;
42class vtkDataSet;
45class vtkFloatArray;
46class vtkGraph;
48class vtkIntArray;
49class vtkPoints;
51
52class VTKIOEXODUS_EXPORT vtkExodusIIReader : public vtkMultiBlockDataSetAlgorithm
53{
54public:
57 void PrintSelf(ostream& os, vtkIndent indent) override;
58
62 virtual int CanReadFile(VTK_FILEPATH const char* fname);
63
64 // virtual void Modified();
65
70
77
79
82 virtual void SetFileName(VTK_FILEPATH const char* fname);
85
87
90 virtual void SetXMLFileName(VTK_FILEPATH const char* fname);
91 vtkGetFilePathMacro(XMLFileName);
93
95
98 vtkSetMacro(TimeStep, int);
99 vtkGetMacro(TimeStep, int);
101
106 void SetModeShape(int val) { this->SetTimeStep(val - 1); }
107
109
115 vtkGetVector2Macro(ModeShapesRange, int);
117
119
124 vtkGetVector2Macro(TimeStepRange, int);
126
128
141 vtkBooleanMacro(GenerateObjectIdCellArray, vtkTypeBool);
142 static const char* GetObjectIdArrayName() { return "ObjectId"; }
144
147 vtkBooleanMacro(GenerateGlobalElementIdArray, vtkTypeBool);
148
151 vtkBooleanMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
152
155 vtkBooleanMacro(GenerateImplicitElementIdArray, vtkTypeBool);
156
159 vtkBooleanMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
160
163 vtkBooleanMacro(GenerateFileIdArray, vtkTypeBool);
164
165 virtual void SetFileId(int f);
167
169
175 enum
176 {
177 SEARCH_TYPE_ELEMENT = 0,
181 ID_NOT_FOUND = -234121312
182 };
183 // NOTE: GetNumberOfObjectTypes must be updated whenever you add an entry to enum ObjectType {...}
185 {
186 // match Exodus macros from exodusII.h and exodusII_ext.h
187 EDGE_BLOCK = 6,
188 FACE_BLOCK = 8,
189 ELEM_BLOCK = 1,
190 NODE_SET = 2,
191 EDGE_SET = 7,
192 FACE_SET = 9,
193 SIDE_SET = 3,
194 ELEM_SET = 10,
195 NODE_MAP = 5,
196 EDGE_MAP = 11,
197 FACE_MAP = 12,
198 ELEM_MAP = 4,
199 GLOBAL = 13,
200 NODAL = 14,
201 // extended values (not in Exodus headers) for use with SetAllArrayStatus:
202 ASSEMBLY = 60,
203 PART = 61,
204 MATERIAL = 62,
205 HIERARCHY = 63,
206 // extended values (not in Exodus headers) for use in cache keys:
207 QA_RECORDS = 103,
208 INFO_RECORDS = 104,
209 GLOBAL_TEMPORAL = 102,
210 NODAL_TEMPORAL = 101,
211 ELEM_BLOCK_TEMPORAL = 100,
212 GLOBAL_CONN = 99,
213 ELEM_BLOCK_ELEM_CONN = 98,
214 ELEM_BLOCK_FACE_CONN =
215 97,
216 ELEM_BLOCK_EDGE_CONN =
217 96,
218 FACE_BLOCK_CONN = 95,
219 EDGE_BLOCK_CONN = 94,
220 ELEM_SET_CONN = 93,
221 SIDE_SET_CONN = 92,
222 FACE_SET_CONN = 91,
223 EDGE_SET_CONN = 90,
224 NODE_SET_CONN = 89,
225 NODAL_COORDS = 88,
226 OBJECT_ID = 87,
227 IMPLICIT_ELEMENT_ID = 108,
228 IMPLICIT_NODE_ID = 107,
229 GLOBAL_ELEMENT_ID =
230 86,
231 GLOBAL_NODE_ID =
232 85,
233 ELEMENT_ID = 84,
234 NODE_ID = 83,
235 NODAL_SQUEEZEMAP = 82,
236 ELEM_BLOCK_ATTRIB = 81,
237 FACE_BLOCK_ATTRIB = 80,
238 EDGE_BLOCK_ATTRIB = 79,
239 FACE_ID = 105,
240 EDGE_ID = 106,
241 ENTITY_COUNTS = 109
242 };
244
245 static const char* GetGlobalElementIdArrayName() { return "GlobalElementId"; }
246 static const char* GetPedigreeElementIdArrayName() { return "PedigreeElementId"; }
247 static int GetGlobalElementID(vtkDataSet* data, int localID);
248 static int GetGlobalElementID(vtkDataSet* data, int localID, int searchType);
249 static const char* GetImplicitElementIdArrayName() { return "ImplicitElementId"; }
250
251 static const char* GetGlobalFaceIdArrayName() { return "GlobalFaceId"; }
252 static const char* GetPedigreeFaceIdArrayName() { return "PedigreeFaceId"; }
253 static int GetGlobalFaceID(vtkDataSet* data, int localID);
254 static int GetGlobalFaceID(vtkDataSet* data, int localID, int searchType);
255 static const char* GetImplicitFaceIdArrayName() { return "ImplicitFaceId"; }
256
257 static const char* GetGlobalEdgeIdArrayName() { return "GlobalEdgeId"; }
258 static const char* GetPedigreeEdgeIdArrayName() { return "PedigreeEdgeId"; }
259 static int GetGlobalEdgeID(vtkDataSet* data, int localID);
260 static int GetGlobalEdgeID(vtkDataSet* data, int localID, int searchType);
261 static const char* GetImplicitEdgeIdArrayName() { return "ImplicitEdgeId"; }
262
264
270 static const char* GetGlobalNodeIdArrayName() { return "GlobalNodeId"; }
271 static const char* GetPedigreeNodeIdArrayName() { return "PedigreeNodeId"; }
272 static int GetGlobalNodeID(vtkDataSet* data, int localID);
273 static int GetGlobalNodeID(vtkDataSet* data, int localID, int searchType);
274 static const char* GetImplicitNodeIdArrayName() { return "ImplicitNodeId"; }
276
281 static const char* GetSideSetSourceElementIdArrayName() { return "SourceElementId"; }
282
287 static const char* GetSideSetSourceElementSideArrayName() { return "SourceElementSide"; }
289
298 vtkBooleanMacro(ApplyDisplacements, vtkTypeBool);
299 virtual void SetDisplacementMagnitude(float s);
302
304
311 vtkBooleanMacro(HasModeShapes, vtkTypeBool);
313
315
322 virtual void SetModeShapeTime(double phase);
325
327
336 vtkBooleanMacro(AnimateModeShapes, vtkTypeBool);
338
340
346 virtual void SetIgnoreFileTime(bool flag);
348 vtkBooleanMacro(IgnoreFileTime, bool);
350
352
355 const char* GetTitle();
359
364
365 int GetObjectTypeFromName(const char* name);
366 const char* GetObjectTypeName(int);
367
369 int GetNumberOfObjects(int objectType);
370 int GetNumberOfEntriesInObject(int objectType, int objectIndex);
371 int GetObjectId(int objectType, int objectIndex);
372 const char* GetObjectName(int objectType, int objectIndex);
373 using Superclass::GetObjectName;
374 int GetObjectIndex(int objectType, const char* objectName);
375 int GetObjectIndex(int objectType, int id);
376 int GetObjectStatus(int objectType, int objectIndex);
377 int GetObjectStatus(int objectType, const char* objectName)
378 {
379 return this->GetObjectStatus(objectType, this->GetObjectIndex(objectType, objectName));
380 }
381 void SetObjectStatus(int objectType, int objectIndex, int status);
382 void SetObjectStatus(int objectType, const char* objectName, int status);
383
385
391 int GetNumberOfObjectArrays(int objectType);
392 const char* GetObjectArrayName(int objectType, int arrayIndex);
393 int GetObjectArrayIndex(int objectType, const char* arrayName);
394 int GetNumberOfObjectArrayComponents(int objectType, int arrayIndex);
395 int GetObjectArrayStatus(int objectType, int arrayIndex);
396 int GetObjectArrayStatus(int objectType, const char* arrayName)
397 {
398 return this->GetObjectArrayStatus(objectType, this->GetObjectArrayIndex(objectType, arrayName));
399 }
400 void SetObjectArrayStatus(int objectType, int arrayIndex, int status);
401 void SetObjectArrayStatus(int objectType, const char* arrayName, int status);
403
405
411 int GetNumberOfObjectAttributes(int objectType, int objectIndex);
412 const char* GetObjectAttributeName(int objectType, int objectIndex, int attribIndex);
413 int GetObjectAttributeIndex(int objectType, int objectIndex, const char* attribName);
414 int GetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex);
415 int GetObjectAttributeStatus(int objectType, int objectIndex, const char* attribName)
416 {
417 return this->GetObjectAttributeStatus(
418 objectType, objectIndex, this->GetObjectAttributeIndex(objectType, objectIndex, attribName));
419 }
420 void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status);
421 void SetObjectAttributeStatus(int objectType, int objectIndex, const char* attribName, int status)
422 {
423 this->SetObjectAttributeStatus(objectType, objectIndex,
424 this->GetObjectAttributeIndex(objectType, objectIndex, attribName), status);
425 }
427
432
434
440 const char* GetPartArrayName(int arrayIdx);
441 int GetPartArrayID(const char* name);
442 const char* GetPartBlockInfo(int arrayIdx);
443 void SetPartArrayStatus(int index, int flag);
444 void SetPartArrayStatus(const char*, int flag);
445 int GetPartArrayStatus(int index);
446 int GetPartArrayStatus(const char*);
448
450
457 const char* GetMaterialArrayName(int arrayIdx);
458 int GetMaterialArrayID(const char* name);
459 void SetMaterialArrayStatus(int index, int flag);
460 void SetMaterialArrayStatus(const char*, int flag);
462 int GetMaterialArrayStatus(const char*);
464
466
473 const char* GetAssemblyArrayName(int arrayIdx);
474 int GetAssemblyArrayID(const char* name);
475 void SetAssemblyArrayStatus(int index, int flag);
476 void SetAssemblyArrayStatus(const char*, int flag);
478 int GetAssemblyArrayStatus(const char*);
480
482
492 const char* GetHierarchyArrayName(int arrayIdx);
493 void SetHierarchyArrayStatus(int index, int flag);
494 void SetHierarchyArrayStatus(const char*, int flag);
496 int GetHierarchyArrayStatus(const char*);
498
499 vtkGetMacro(DisplayType, int);
500 virtual void SetDisplayType(int type);
501
505 int IsValidVariable(const char* type, const char* name);
506
510 int GetVariableID(const char* type, const char* name);
511
512 void SetAllArrayStatus(int otype, int status);
513 // Helper functions
514 // static int StringsEqual(const char* s1, char* s2);
515 // static void StringUppercase(const char* str, char* upperstr);
516 // static char *StrDupWithNew(const char *s);
517
518 // time series query functions
519 int GetTimeSeriesData(int ID, const char* vName, const char* vType, vtkFloatArray* result);
520
521 int GetNumberOfEdgeBlockArrays() { return this->GetNumberOfObjects(EDGE_BLOCK); }
522 const char* GetEdgeBlockArrayName(int index) { return this->GetObjectName(EDGE_BLOCK, index); }
523 int GetEdgeBlockArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_BLOCK, name); }
524 void SetEdgeBlockArrayStatus(const char* name, int flag)
525 {
526 this->SetObjectStatus(EDGE_BLOCK, name, flag);
527 }
528
529 int GetNumberOfFaceBlockArrays() { return this->GetNumberOfObjects(FACE_BLOCK); }
530 const char* GetFaceBlockArrayName(int index) { return this->GetObjectName(FACE_BLOCK, index); }
531 int GetFaceBlockArrayStatus(const char* name) { return this->GetObjectStatus(FACE_BLOCK, name); }
532 void SetFaceBlockArrayStatus(const char* name, int flag)
533 {
534 this->SetObjectStatus(FACE_BLOCK, name, flag);
535 }
536
537 int GetNumberOfElementBlockArrays() { return this->GetNumberOfObjects(ELEM_BLOCK); }
538 const char* GetElementBlockArrayName(int index) { return this->GetObjectName(ELEM_BLOCK, index); }
539 int GetElementBlockArrayStatus(const char* name)
540 {
541 return this->GetObjectStatus(ELEM_BLOCK, name);
542 }
543 void SetElementBlockArrayStatus(const char* name, int flag)
544 {
545 this->SetObjectStatus(ELEM_BLOCK, name, flag);
546 }
547
548 int GetNumberOfGlobalResultArrays() { return this->GetNumberOfObjectArrays(GLOBAL); }
549 const char* GetGlobalResultArrayName(int index)
550 {
551 return this->GetObjectArrayName(GLOBAL, index);
552 }
553 int GetGlobalResultArrayStatus(const char* name)
554 {
555 return this->GetObjectArrayStatus(GLOBAL, name);
556 }
557 void SetGlobalResultArrayStatus(const char* name, int flag)
558 {
559 this->SetObjectArrayStatus(GLOBAL, name, flag);
560 }
561
562 int GetNumberOfPointResultArrays() { return this->GetNumberOfObjectArrays(NODAL); }
563 const char* GetPointResultArrayName(int index) { return this->GetObjectArrayName(NODAL, index); }
564 int GetPointResultArrayStatus(const char* name)
565 {
566 return this->GetObjectArrayStatus(NODAL, name);
567 }
568 void SetPointResultArrayStatus(const char* name, int flag)
569 {
570 this->SetObjectArrayStatus(NODAL, name, flag);
571 }
572
573 int GetNumberOfEdgeResultArrays() { return this->GetNumberOfObjectArrays(EDGE_BLOCK); }
574 const char* GetEdgeResultArrayName(int index)
575 {
576 return this->GetObjectArrayName(EDGE_BLOCK, index);
577 }
578 int GetEdgeResultArrayStatus(const char* name)
579 {
580 return this->GetObjectArrayStatus(EDGE_BLOCK, name);
581 }
582 void SetEdgeResultArrayStatus(const char* name, int flag)
583 {
584 this->SetObjectArrayStatus(EDGE_BLOCK, name, flag);
585 }
586
587 int GetNumberOfFaceResultArrays() { return this->GetNumberOfObjectArrays(FACE_BLOCK); }
588 const char* GetFaceResultArrayName(int index)
589 {
590 return this->GetObjectArrayName(FACE_BLOCK, index);
591 }
592 int GetFaceResultArrayStatus(const char* name)
593 {
594 return this->GetObjectArrayStatus(FACE_BLOCK, name);
595 }
596 void SetFaceResultArrayStatus(const char* name, int flag)
597 {
598 this->SetObjectArrayStatus(FACE_BLOCK, name, flag);
599 }
600
601 int GetNumberOfElementResultArrays() { return this->GetNumberOfObjectArrays(ELEM_BLOCK); }
602 const char* GetElementResultArrayName(int index)
603 {
604 return this->GetObjectArrayName(ELEM_BLOCK, index);
605 }
606 int GetElementResultArrayStatus(const char* name)
607 {
608 return this->GetObjectArrayStatus(ELEM_BLOCK, name);
609 }
610 void SetElementResultArrayStatus(const char* name, int flag)
611 {
612 this->SetObjectArrayStatus(ELEM_BLOCK, name, flag);
613 }
614
615 int GetNumberOfNodeMapArrays() { return this->GetNumberOfObjects(NODE_MAP); }
616 const char* GetNodeMapArrayName(int index) { return this->GetObjectName(NODE_MAP, index); }
617 int GetNodeMapArrayStatus(const char* name) { return this->GetObjectStatus(NODE_MAP, name); }
618 void SetNodeMapArrayStatus(const char* name, int flag)
619 {
620 this->SetObjectStatus(NODE_MAP, name, flag);
621 }
622
623 int GetNumberOfEdgeMapArrays() { return this->GetNumberOfObjects(EDGE_MAP); }
624 const char* GetEdgeMapArrayName(int index) { return this->GetObjectName(EDGE_MAP, index); }
625 int GetEdgeMapArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_MAP, name); }
626 void SetEdgeMapArrayStatus(const char* name, int flag)
627 {
628 this->SetObjectStatus(EDGE_MAP, name, flag);
629 }
630
631 int GetNumberOfFaceMapArrays() { return this->GetNumberOfObjects(FACE_MAP); }
632 const char* GetFaceMapArrayName(int index) { return this->GetObjectName(FACE_MAP, index); }
633 int GetFaceMapArrayStatus(const char* name) { return this->GetObjectStatus(FACE_MAP, name); }
634 void SetFaceMapArrayStatus(const char* name, int flag)
635 {
636 this->SetObjectStatus(FACE_MAP, name, flag);
637 }
638
639 int GetNumberOfElementMapArrays() { return this->GetNumberOfObjects(ELEM_MAP); }
640 const char* GetElementMapArrayName(int index) { return this->GetObjectName(ELEM_MAP, index); }
641 int GetElementMapArrayStatus(const char* name) { return this->GetObjectStatus(ELEM_MAP, name); }
642 void SetElementMapArrayStatus(const char* name, int flag)
643 {
644 this->SetObjectStatus(ELEM_MAP, name, flag);
645 }
646
647 int GetNumberOfNodeSetArrays() { return this->GetNumberOfObjects(NODE_SET); }
648 const char* GetNodeSetArrayName(int index) { return this->GetObjectName(NODE_SET, index); }
649 int GetNodeSetArrayStatus(const char* name) { return this->GetObjectStatus(NODE_SET, name); }
650 void SetNodeSetArrayStatus(const char* name, int flag)
651 {
652 this->SetObjectStatus(NODE_SET, name, flag);
653 }
654
655 int GetNumberOfSideSetArrays() { return this->GetNumberOfObjects(SIDE_SET); }
656 const char* GetSideSetArrayName(int index) { return this->GetObjectName(SIDE_SET, index); }
657 int GetSideSetArrayStatus(const char* name) { return this->GetObjectStatus(SIDE_SET, name); }
658 void SetSideSetArrayStatus(const char* name, int flag)
659 {
660 this->SetObjectStatus(SIDE_SET, name, flag);
661 }
662
663 int GetNumberOfEdgeSetArrays() { return this->GetNumberOfObjects(EDGE_SET); }
664 const char* GetEdgeSetArrayName(int index) { return this->GetObjectName(EDGE_SET, index); }
665 int GetEdgeSetArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_SET, name); }
666 void SetEdgeSetArrayStatus(const char* name, int flag)
667 {
668 this->SetObjectStatus(EDGE_SET, name, flag);
669 }
670
671 int GetNumberOfFaceSetArrays() { return this->GetNumberOfObjects(FACE_SET); }
672 const char* GetFaceSetArrayName(int index) { return this->GetObjectName(FACE_SET, index); }
673 int GetFaceSetArrayStatus(const char* name) { return this->GetObjectStatus(FACE_SET, name); }
674 void SetFaceSetArrayStatus(const char* name, int flag)
675 {
676 this->SetObjectStatus(FACE_SET, name, flag);
677 }
678
679 int GetNumberOfElementSetArrays() { return this->GetNumberOfObjects(ELEM_SET); }
680 const char* GetElementSetArrayName(int index) { return this->GetObjectName(ELEM_SET, index); }
681 int GetElementSetArrayStatus(const char* name) { return this->GetObjectStatus(ELEM_SET, name); }
682 void SetElementSetArrayStatus(const char* name, int flag)
683 {
684 this->SetObjectStatus(ELEM_SET, name, flag);
685 }
686
687 int GetNumberOfNodeSetResultArrays() { return this->GetNumberOfObjectArrays(NODE_SET); }
688 const char* GetNodeSetResultArrayName(int index)
689 {
690 return this->GetObjectArrayName(NODE_SET, index);
691 }
692 int GetNodeSetResultArrayStatus(const char* name)
693 {
694 return this->GetObjectArrayStatus(NODE_SET, name);
695 }
696 void SetNodeSetResultArrayStatus(const char* name, int flag)
697 {
698 this->SetObjectArrayStatus(NODE_SET, name, flag);
699 }
700
701 int GetNumberOfSideSetResultArrays() { return this->GetNumberOfObjectArrays(SIDE_SET); }
702 const char* GetSideSetResultArrayName(int index)
703 {
704 return this->GetObjectArrayName(SIDE_SET, index);
705 }
706 int GetSideSetResultArrayStatus(const char* name)
707 {
708 return this->GetObjectArrayStatus(SIDE_SET, name);
709 }
710 void SetSideSetResultArrayStatus(const char* name, int flag)
711 {
712 this->SetObjectArrayStatus(SIDE_SET, name, flag);
713 }
714
715 int GetNumberOfEdgeSetResultArrays() { return this->GetNumberOfObjectArrays(EDGE_SET); }
716 const char* GetEdgeSetResultArrayName(int index)
717 {
718 return this->GetObjectArrayName(EDGE_SET, index);
719 }
720 int GetEdgeSetResultArrayStatus(const char* name)
721 {
722 return this->GetObjectArrayStatus(EDGE_SET, name);
723 }
724 void SetEdgeSetResultArrayStatus(const char* name, int flag)
725 {
726 this->SetObjectArrayStatus(EDGE_SET, name, flag);
727 }
728
729 int GetNumberOfFaceSetResultArrays() { return this->GetNumberOfObjectArrays(FACE_SET); }
730 const char* GetFaceSetResultArrayName(int index)
731 {
732 return this->GetObjectArrayName(FACE_SET, index);
733 }
734 int GetFaceSetResultArrayStatus(const char* name)
735 {
736 return this->GetObjectArrayStatus(FACE_SET, name);
737 }
738 void SetFaceSetResultArrayStatus(const char* name, int flag)
739 {
740 this->SetObjectArrayStatus(FACE_SET, name, flag);
741 }
742
743 int GetNumberOfElementSetResultArrays() { return this->GetNumberOfObjectArrays(ELEM_SET); }
744 const char* GetElementSetResultArrayName(int index)
745 {
746 return this->GetObjectArrayName(ELEM_SET, index);
747 }
748 int GetElementSetResultArrayStatus(const char* name)
749 {
750 return this->GetObjectArrayStatus(ELEM_SET, name);
751 }
752 void SetElementSetResultArrayStatus(const char* name, int flag)
753 {
754 this->SetObjectArrayStatus(ELEM_SET, name, flag);
755 }
756
765 void Reset();
766
776
781
785 void SetCacheSize(double CacheSize);
786
790 double GetCacheSize();
791
793
805 void SetSqueezePoints(bool sp);
808
809 virtual void Dump();
810
816
818
821 vtkGetMacro(SILUpdateStamp, int);
823
825
831
833
844
846
853 vtkSetMacro(UseLegacyBlockNames, bool);
854 vtkGetMacro(UseLegacyBlockNames, bool);
855 vtkBooleanMacro(UseLegacyBlockNames, bool);
857protected:
860
861 // helper for finding IDs
862 static int GetIDHelper(const char* arrayName, vtkDataSet* data, int localID, int searchType);
863 static int GetGlobalID(const char* arrayName, vtkDataSet* data, int localID, int searchType);
864
866 vtkGetObjectMacro(Metadata, vtkExodusIIReaderPrivate);
867
873
874 // Time query function. Called by ExecuteInformation().
875 // Fills the TimestepValues array.
877
882
887 // int RequestDataOverTime( vtkInformation *, vtkInformationVector **, vtkInformationVector *);
888
889 // Parameters for controlling what is read in.
890 char* FileName;
893 int TimeStepRange[2];
896
897 // Information specific for exodus files.
898
899 // 1=display Block names
900 // 2=display Part names
901 // 3=display Material names
903
904 // Metadata containing a description of the currently open file.
906
908
909 friend class vtkPExodusIIReader;
910
911private:
912 vtkExodusIIReader(const vtkExodusIIReader&) = delete;
913 void operator=(const vtkExodusIIReader&) = delete;
914
915 void AddDisplacements(vtkUnstructuredGrid* output);
916 int ModeShapesRange[2];
917
918 bool UseLegacyBlockNames;
919};
920
921VTK_ABI_NAMESPACE_END
922#endif
abstract superclass for arrays of numeric data
abstract class to specify dataset behavior
Definition vtkDataSet.h:62
This class holds metadata for an Exodus file.
Read exodus 2 files .ex2.
virtual void SetGenerateGlobalNodeIdArray(vtkTypeBool g)
static const char * GetSideSetSourceElementIdArrayName()
Get the name of the array that stores the mapping from side set cells back to the global id of the el...
int GetNumberOfElementsInFile()
int IsValidVariable(const char *type, const char *name)
return boolean indicating whether the type,name is a valid variable
vtkTypeBool ProcessRequest(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
Upstream/Downstream requests form the generalized interface through which executives invoke a algorit...
vtkTypeBool GetAnimateModeShapes()
If this flag is on (the default) and HasModeShapes is also on, then this reader will report a continu...
virtual void SetFileName(VTK_FILEPATH const char *fname)
Specify file name of the Exodus file.
const char * GetObjectTypeName(int)
const char * GetNodeSetArrayName(int index)
int GetEdgeBlockArrayStatus(const char *name)
int GetFaceResultArrayStatus(const char *name)
virtual void SetFileId(int f)
void SetEdgeBlockArrayStatus(const char *name, int flag)
int GetNumberOfFacesInFile()
static int GetGlobalNodeID(vtkDataSet *data, int localID, int searchType)
Extra point data array that can be generated.
vtkTypeBool GetGenerateImplicitNodeIdArray()
static int GetGlobalFaceID(vtkDataSet *data, int localID)
int GetObjectArrayStatus(int objectType, int arrayIndex)
By default arrays are not loaded.
static const char * GetImplicitNodeIdArrayName()
Extra point data array that can be generated.
int GetObjectIndex(int objectType, const char *objectName)
void SetElementResultArrayStatus(const char *name, int flag)
int GetNumberOfObjectArrays(int objectType)
By default arrays are not loaded.
int GetEdgeSetResultArrayStatus(const char *name)
static const char * GetImplicitFaceIdArrayName()
void SetElementMapArrayStatus(const char *name, int flag)
void SetElementSetArrayStatus(const char *name, int flag)
const char * GetFaceResultArrayName(int index)
void SetSideSetResultArrayStatus(const char *name, int flag)
int GetMaterialArrayStatus(const char *)
By default all materials are loaded.
static vtkInformationIntegerKey * GLOBAL_VARIABLE()
Exodus reader outputs global variables and global temporal variables, together with some other variab...
int GetNodeMapArrayStatus(const char *name)
int GetNumberOfHierarchyArrays()
By default all hierarchy entries are loaded.
static int GetGlobalEdgeID(vtkDataSet *data, int localID, int searchType)
int GetElementMapArrayStatus(const char *name)
void SetEdgeResultArrayStatus(const char *name, int flag)
static const char * GetGlobalEdgeIdArrayName()
int GetNumberOfPartArrays()
By default all parts are loaded.
int GetHierarchyArrayStatus(const char *)
By default all hierarchy entries are loaded.
int GetNumberOfTimeSteps()
Access to meta data generated by UpdateInformation.
vtkTypeBool GetGenerateFileIdArray()
int GetNumberOfNodesInFile()
virtual void Dump()
const char * GetEdgeBlockArrayName(int index)
void SetFaceBlockArrayStatus(const char *name, int flag)
int GetPointResultArrayStatus(const char *name)
int GetPartArrayStatus(int index)
By default all parts are loaded.
const char * GetFaceMapArrayName(int index)
void SetPartArrayStatus(int index, int flag)
By default all parts are loaded.
virtual void SetHasModeShapes(vtkTypeBool ms)
Set/Get whether the Exodus sequence number corresponds to time steps or mode shapes.
static const char * GetPedigreeFaceIdArrayName()
static int GetGlobalID(const char *arrayName, vtkDataSet *data, int localID, int searchType)
vtkTypeBool GetGenerateGlobalElementIdArray()
virtual void SetGenerateImplicitNodeIdArray(vtkTypeBool g)
const char * GetSideSetArrayName(int index)
int GetPartArrayID(const char *name)
By default all parts are loaded.
~vtkExodusIIReader() override
bool GetSqueezePoints()
Should the reader output only points used by elements in the output mesh, or all the points.
const char * GetObjectAttributeName(int objectType, int objectIndex, int attribIndex)
By default attributes are not loaded.
int GetEdgeMapArrayStatus(const char *name)
void SetPartArrayStatus(const char *, int flag)
By default all parts are loaded.
vtkExodusIIReaderPrivate * Metadata
int GetMaterialArrayStatus(int index)
By default all materials are loaded.
int GetElementResultArrayStatus(const char *name)
int GetNumberOfAssemblyArrays()
By default all assemblies are loaded.
void SetAssemblyArrayStatus(const char *, int flag)
By default all assemblies are loaded.
int GetEdgeSetArrayStatus(const char *name)
void SetNodeSetResultArrayStatus(const char *name, int flag)
void SetMaterialArrayStatus(int index, int flag)
By default all materials are loaded.
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
const char * GetFaceSetArrayName(int index)
const char * GetGlobalResultArrayName(int index)
int GetNodeSetResultArrayStatus(const char *name)
virtual void SetGenerateGlobalElementIdArray(vtkTypeBool g)
void SetAllArrayStatus(int otype, int status)
int GetHierarchyArrayStatus(int index)
By default all hierarchy entries are loaded.
vtkTypeBool GetGenerateObjectIdCellArray()
Extra cell data array that can be generated.
int GetNumberOfMaterialArrays()
By default all materials are loaded.
vtkTypeBool GetGenerateImplicitElementIdArray()
const char * GetEdgeMapArrayName(int index)
void ResetCache()
Clears out the cache entries.
virtual void SetMetadata(vtkExodusIIReaderPrivate *)
void SetMaterialArrayStatus(const char *, int flag)
By default all materials are loaded.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
static const char * GetGlobalFaceIdArrayName()
void Reset()
Reset the user-specified parameters and flush internal arrays so that the reader state is just as it ...
double GetCacheSize()
Get the size of the cache in MiB.
virtual void SetDisplayType(int type)
int GetNumberOfEntriesInObject(int objectType, int objectIndex)
int GetObjectStatus(int objectType, int objectIndex)
void SetSideSetArrayStatus(const char *name, int flag)
void SetNodeMapArrayStatus(const char *name, int flag)
static int GetIDHelper(const char *arrayName, vtkDataSet *data, int localID, int searchType)
void SetEdgeMapArrayStatus(const char *name, int flag)
void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status)
By default attributes are not loaded.
vtkTimeStamp FileNameMTime
int GetAssemblyArrayStatus(const char *)
By default all assemblies are loaded.
const char * GetFaceSetResultArrayName(int index)
void SetModeShape(int val)
Convenience method to set the mode-shape which is same as this->SetTimeStep(val-1);.
static const char * GetPedigreeEdgeIdArrayName()
int GetAssemblyArrayID(const char *name)
By default all assemblies are loaded.
static int GetGlobalEdgeID(vtkDataSet *data, int localID)
void SetObjectArrayStatus(int objectType, int arrayIndex, int status)
By default arrays are not loaded.
void SetFaceMapArrayStatus(const char *name, int flag)
bool GetIgnoreFileTime()
When on, this option ignores the time values assigned to each time step in the file.
int GetSideSetArrayStatus(const char *name)
void ResetSettings()
Reset the user-specified parameters to their default values.
int GetMaterialArrayID(const char *name)
By default all materials are loaded.
const char * GetHierarchyArrayName(int arrayIdx)
By default all hierarchy entries are loaded.
void SetFaceResultArrayStatus(const char *name, int flag)
int GetObjectIndex(int objectType, int id)
int GetObjectAttributeStatus(int objectType, int objectIndex, const char *attribName)
By default attributes are not loaded.
static int GetGlobalFaceID(vtkDataSet *data, int localID, int searchType)
void SetObjectStatus(int objectType, int objectIndex, int status)
const char * GetElementMapArrayName(int index)
void SetElementSetResultArrayStatus(const char *name, int flag)
int GetDimensionality()
Access to meta data generated by UpdateInformation.
void SetAssemblyArrayStatus(int index, int flag)
By default all assemblies are loaded.
int GetObjectArrayIndex(int objectType, const char *arrayName)
By default arrays are not loaded.
static const char * GetPedigreeElementIdArrayName()
int GetNumberOfObjectAttributes(int objectType, int objectIndex)
By default attributes are not loaded.
void GetAllTimes(vtkInformationVector *)
int GetFaceSetArrayStatus(const char *name)
void SetHierarchyArrayStatus(int index, int flag)
By default all hierarchy entries are loaded.
int GetMaxNameLength()
Get the max_name_length in the file.
static const char * GetPedigreeNodeIdArrayName()
Extra point data array that can be generated.
static const char * GetImplicitElementIdArrayName()
double GetModeShapeTime()
Set/Get the time used to animate mode shapes.
vtkTimeStamp XMLFileNameMTime
int GetEdgeResultArrayStatus(const char *name)
float GetDisplacementMagnitude()
Geometric locations can include displacements.
virtual void SetGenerateObjectIdCellArray(vtkTypeBool g)
Extra cell data array that can be generated.
int GetTimeSeriesData(int ID, const char *vName, const char *vType, vtkFloatArray *result)
virtual int CanReadFile(VTK_FILEPATH const char *fname)
Determine if the file can be read with this reader.
void SetObjectStatus(int objectType, const char *objectName, int status)
void SetObjectAttributeStatus(int objectType, int objectIndex, const char *attribName, int status)
By default attributes are not loaded.
static int GetGlobalNodeID(vtkDataSet *data, int localID)
Extra point data array that can be generated.
int GetFaceMapArrayStatus(const char *name)
int GetFaceBlockArrayStatus(const char *name)
const char * GetSideSetResultArrayName(int index)
virtual vtkIdType GetTotalNumberOfEdges()
virtual vtkMTimeType GetMetadataMTime()
Return the MTime of the internal data structure.
virtual void SetXMLFileName(VTK_FILEPATH const char *fname)
Specify file name of the xml file.
bool FindXMLFile()
Returns true if the file given by XMLFileName exists.
void AdvertiseTimeSteps(vtkInformation *outputInfo)
Populates the TIME_STEPS and TIME_RANGE keys based on file metadata.
const char * GetNodeMapArrayName(int index)
const char * GetMaterialArrayName(int arrayIdx)
By default all materials are loaded.
vtkTypeBool GetHasModeShapes()
Set/Get whether the Exodus sequence number corresponds to time steps or mode shapes.
void SetEdgeSetArrayStatus(const char *name, int flag)
virtual vtkIdType GetTotalNumberOfFaces()
vtkGraph * GetSIL()
SIL describes organization of/relationships between classifications eg.
int GetObjectId(int objectType, int objectIndex)
int GetNumberOfElementSetResultArrays()
const char * GetTitle()
Access to meta data generated by UpdateInformation.
int GetNumberOfEdgesInFile()
void SetElementBlockArrayStatus(const char *name, int flag)
static const char * GetGlobalElementIdArrayName()
virtual void SetModeShapeTime(double phase)
Set/Get the time used to animate mode shapes.
const char * GetPartBlockInfo(int arrayIdx)
By default all parts are loaded.
void SetFaceSetArrayStatus(const char *name, int flag)
void SetGlobalResultArrayStatus(const char *name, int flag)
vtkTypeBool GetGenerateGlobalNodeIdArray()
const char * GetPartArrayName(int arrayIdx)
By default all parts are loaded.
vtkGetFilePathMacro(XMLFileName)
Specify file name of the xml file.
virtual void SetApplyDisplacements(vtkTypeBool d)
Geometric locations can include displacements.
static const char * GetSideSetSourceElementSideArrayName()
Get the name of the array that stores the mapping from side set cells back to the canonical side of t...
virtual vtkIdType GetTotalNumberOfElements()
const char * GetNodeSetResultArrayName(int index)
vtkMTimeType GetMTime() override
Return the object's MTime.
void SetHierarchyArrayStatus(const char *, int flag)
By default all hierarchy entries are loaded.
static const char * GetGlobalNodeIdArrayName()
Extra point data array that can be generated.
int GetNumberOfObjects(int objectType)
virtual void SetIgnoreFileTime(bool flag)
When on, this option ignores the time values assigned to each time step in the file.
int GetAssemblyArrayStatus(int index)
By default all assemblies are loaded.
int GetObjectTypeFromName(const char *name)
int GetElementBlockArrayStatus(const char *name)
const char * GetEdgeResultArrayName(int index)
int GetElementSetResultArrayStatus(const char *name)
int GetObjectAttributeIndex(int objectType, int objectIndex, const char *attribName)
By default attributes are not loaded.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
const char * GetElementResultArrayName(int index)
virtual void SetDisplacementMagnitude(float s)
Geometric locations can include displacements.
int GetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex)
By default attributes are not loaded.
int GetGlobalResultArrayStatus(const char *name)
virtual vtkIdType GetTotalNumberOfNodes()
virtual void SetAnimateModeShapes(vtkTypeBool flag)
If this flag is on (the default) and HasModeShapes is also on, then this reader will report a continu...
virtual void SetGenerateImplicitElementIdArray(vtkTypeBool g)
const char * GetEdgeSetArrayName(int index)
void SetPointResultArrayStatus(const char *name, int flag)
const char * GetFaceBlockArrayName(int index)
static const char * GetImplicitEdgeIdArrayName()
int GetVariableID(const char *type, const char *name)
Return the id of the type,name variable.
void SetNodeSetArrayStatus(const char *name, int flag)
int GetSideSetResultArrayStatus(const char *name)
void SetSqueezePoints(bool sp)
Should the reader output only points used by elements in the output mesh, or all the points.
virtual void SetGenerateFileIdArray(vtkTypeBool f)
void SetFaceSetResultArrayStatus(const char *name, int flag)
const char * GetObjectName(int objectType, int objectIndex)
void SetCacheSize(double CacheSize)
Set the size of the cache in MiB.
static int GetGlobalElementID(vtkDataSet *data, int localID, int searchType)
int GetElementSetArrayStatus(const char *name)
vtkGetFilePathMacro(FileName)
Specify file name of the Exodus file.
const char * GetObjectArrayName(int objectType, int arrayIndex)
By default arrays are not loaded.
ObjectType
Extra cell data array that can be generated.
int GetNumberOfObjectArrayComponents(int objectType, int arrayIndex)
By default arrays are not loaded.
static const char * GetObjectIdArrayName()
Extra cell data array that can be generated.
vtkTypeBool GetApplyDisplacements()
Geometric locations can include displacements.
int GetObjectStatus(int objectType, const char *objectName)
const char * GetPointResultArrayName(int index)
int GetObjectArrayStatus(int objectType, const char *arrayName)
By default arrays are not loaded.
void SetEdgeSetResultArrayStatus(const char *name, int flag)
static int GetGlobalElementID(vtkDataSet *data, int localID)
static vtkExodusIIReader * New()
int GetNodeSetArrayStatus(const char *name)
const char * GetElementSetResultArrayName(int index)
const char * GetEdgeSetResultArrayName(int index)
int GetFaceSetResultArrayStatus(const char *name)
const char * GetAssemblyArrayName(int arrayIdx)
By default all assemblies are loaded.
const char * GetElementSetArrayName(int index)
int GetPartArrayStatus(const char *)
By default all parts are loaded.
static vtkInformationIntegerKey * GLOBAL_TEMPORAL_VARIABLE()
Exodus reader outputs global variables and global temporal variables, together with some other variab...
void SetObjectArrayStatus(int objectType, const char *arrayName, int status)
By default arrays are not loaded.
const char * GetElementBlockArrayName(int index)
dynamic, self-adjusting array of float
Base class for graph data types.
Definition vtkGraph.h:290
a simple class to control print indentation
Definition vtkIndent.h:38
Key for integer values in vtkInformation.
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
Definition vtkIntArray.h:44
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
virtual std::string GetObjectName() const
Set/get the name of this object for reporting purposes.
Read Exodus II files (.exii)
represent and manipulate 3D points
Definition vtkPoints.h:38
record modification and/or execution time
dataset represents arbitrary combinations of all possible cell types
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270
#define VTK_FILEPATH