Thanks to visit codestin.com
Credit goes to vtk.org

VTK  9.5.20250922
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
54#ifndef vtkExodusIIReader_h
55#define vtkExodusIIReader_h
56
57#include "vtkIOExodusModule.h" // For export macro
59
60VTK_ABI_NAMESPACE_BEGIN
61class vtkDataArray;
62class vtkDataSet;
65class vtkFloatArray;
66class vtkGraph;
68class vtkIntArray;
69class vtkPoints;
71
72class VTKIOEXODUS_EXPORT vtkExodusIIReader : public vtkMultiBlockDataSetAlgorithm
73{
74public:
77 void PrintSelf(ostream& os, vtkIndent indent) override;
78
82 virtual int CanReadFile(VTK_FILEPATH const char* fname);
83
84 // virtual void Modified();
85
90
97
99
102 virtual void SetFileName(VTK_FILEPATH const char* fname);
105
107
110 virtual void SetXMLFileName(VTK_FILEPATH const char* fname);
113
115
118 vtkSetMacro(TimeStep, int);
119 vtkGetMacro(TimeStep, int);
121
126 void SetModeShape(int val) { this->SetTimeStep(val - 1); }
127
129
135 vtkGetVector2Macro(ModeShapesRange, int);
137
139
144 vtkGetVector2Macro(TimeStepRange, int);
146
148
161 vtkBooleanMacro(GenerateObjectIdCellArray, vtkTypeBool);
162 static const char* GetObjectIdArrayName() { return "ObjectId"; }
164
167 vtkBooleanMacro(GenerateGlobalElementIdArray, vtkTypeBool);
168
171 vtkBooleanMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
172
175 vtkBooleanMacro(GenerateImplicitElementIdArray, vtkTypeBool);
176
179 vtkBooleanMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
180
183 vtkBooleanMacro(GenerateFileIdArray, vtkTypeBool);
184
185 virtual void SetFileId(int f);
187
189
195 // NOLINTNEXTLINE(readability-enum-initial-value)
196 enum
197 {
198 SEARCH_TYPE_ELEMENT = 0,
202 ID_NOT_FOUND = -234121312
203 };
204 // NOTE: GetNumberOfObjectTypes must be updated whenever you add an entry to enum ObjectType {...}
206 {
207 // match Exodus macros from exodusII.h and exodusII_ext.h
208 EDGE_BLOCK = 6,
209 FACE_BLOCK = 8,
210 ELEM_BLOCK = 1,
211 NODE_SET = 2,
212 EDGE_SET = 7,
213 FACE_SET = 9,
214 SIDE_SET = 3,
215 ELEM_SET = 10,
216 NODE_MAP = 5,
217 EDGE_MAP = 11,
218 FACE_MAP = 12,
219 ELEM_MAP = 4,
220 GLOBAL = 13,
221 NODAL = 14,
222 // extended values (not in Exodus headers) for use with SetAllArrayStatus:
223 ASSEMBLY = 60,
224 PART = 61,
225 MATERIAL = 62,
226 HIERARCHY = 63,
227 // extended values (not in Exodus headers) for use in cache keys:
228 QA_RECORDS = 103,
229 INFO_RECORDS = 104,
230 GLOBAL_TEMPORAL = 102,
231 NODAL_TEMPORAL = 101,
232 ELEM_BLOCK_TEMPORAL = 100,
233 GLOBAL_CONN = 99,
234 ELEM_BLOCK_ELEM_CONN = 98,
235 ELEM_BLOCK_FACE_CONN =
236 97,
237 ELEM_BLOCK_EDGE_CONN =
238 96,
239 FACE_BLOCK_CONN = 95,
240 EDGE_BLOCK_CONN = 94,
241 ELEM_SET_CONN = 93,
242 SIDE_SET_CONN = 92,
243 FACE_SET_CONN = 91,
244 EDGE_SET_CONN = 90,
245 NODE_SET_CONN = 89,
246 NODAL_COORDS = 88,
247 OBJECT_ID = 87,
248 IMPLICIT_ELEMENT_ID = 108,
249 IMPLICIT_NODE_ID = 107,
250 GLOBAL_ELEMENT_ID =
251 86,
252 GLOBAL_NODE_ID =
253 85,
254 ELEMENT_ID = 84,
255 NODE_ID = 83,
256 NODAL_SQUEEZEMAP = 82,
257 ELEM_BLOCK_ATTRIB = 81,
258 FACE_BLOCK_ATTRIB = 80,
259 EDGE_BLOCK_ATTRIB = 79,
260 FACE_ID = 105,
261 EDGE_ID = 106,
262 ENTITY_COUNTS = 109
263 };
265
266 static const char* GetGlobalElementIdArrayName() { return "GlobalElementId"; }
267 static const char* GetPedigreeElementIdArrayName() { return "PedigreeElementId"; }
268 static int GetGlobalElementID(vtkDataSet* data, int localID);
269 static int GetGlobalElementID(vtkDataSet* data, int localID, int searchType);
270 static const char* GetImplicitElementIdArrayName() { return "ImplicitElementId"; }
271
272 static const char* GetGlobalFaceIdArrayName() { return "GlobalFaceId"; }
273 static const char* GetPedigreeFaceIdArrayName() { return "PedigreeFaceId"; }
274 static int GetGlobalFaceID(vtkDataSet* data, int localID);
275 static int GetGlobalFaceID(vtkDataSet* data, int localID, int searchType);
276 static const char* GetImplicitFaceIdArrayName() { return "ImplicitFaceId"; }
277
278 static const char* GetGlobalEdgeIdArrayName() { return "GlobalEdgeId"; }
279 static const char* GetPedigreeEdgeIdArrayName() { return "PedigreeEdgeId"; }
280 static int GetGlobalEdgeID(vtkDataSet* data, int localID);
281 static int GetGlobalEdgeID(vtkDataSet* data, int localID, int searchType);
282 static const char* GetImplicitEdgeIdArrayName() { return "ImplicitEdgeId"; }
283
285
291 static const char* GetGlobalNodeIdArrayName() { return "GlobalNodeId"; }
292 static const char* GetPedigreeNodeIdArrayName() { return "PedigreeNodeId"; }
293 static int GetGlobalNodeID(vtkDataSet* data, int localID);
294 static int GetGlobalNodeID(vtkDataSet* data, int localID, int searchType);
295 static const char* GetImplicitNodeIdArrayName() { return "ImplicitNodeId"; }
297
302 static const char* GetSideSetSourceElementIdArrayName() { return "SourceElementId"; }
303
308 static const char* GetSideSetSourceElementSideArrayName() { return "SourceElementSide"; }
310
319 vtkBooleanMacro(ApplyDisplacements, vtkTypeBool);
320 virtual void SetDisplacementMagnitude(float s);
323
325
332 vtkBooleanMacro(HasModeShapes, vtkTypeBool);
334
336
343 virtual void SetModeShapeTime(double phase);
346
348
357 vtkBooleanMacro(AnimateModeShapes, vtkTypeBool);
359
361
367 virtual void SetIgnoreFileTime(bool flag);
369 vtkBooleanMacro(IgnoreFileTime, bool);
371
373
376 const char* GetTitle();
380
385
386 int GetObjectTypeFromName(const char* name);
387 const char* GetObjectTypeName(int);
388
390 int GetNumberOfObjects(int objectType);
391 int GetNumberOfEntriesInObject(int objectType, int objectIndex);
392 int GetObjectId(int objectType, int objectIndex);
393 const char* GetObjectName(int objectType, int objectIndex);
394 using Superclass::GetObjectName;
395 int GetObjectIndex(int objectType, const char* objectName);
396 int GetObjectIndex(int objectType, int id);
397 int GetObjectStatus(int objectType, int objectIndex);
398 int GetObjectStatus(int objectType, const char* objectName)
399 {
400 return this->GetObjectStatus(objectType, this->GetObjectIndex(objectType, objectName));
401 }
402 void SetObjectStatus(int objectType, int objectIndex, int status);
403 void SetObjectStatus(int objectType, const char* objectName, int status);
404
406
412 int GetNumberOfObjectArrays(int objectType);
413 const char* GetObjectArrayName(int objectType, int arrayIndex);
414 int GetObjectArrayIndex(int objectType, const char* arrayName);
415 int GetNumberOfObjectArrayComponents(int objectType, int arrayIndex);
416 int GetObjectArrayStatus(int objectType, int arrayIndex);
417 int GetObjectArrayStatus(int objectType, const char* arrayName)
418 {
419 return this->GetObjectArrayStatus(objectType, this->GetObjectArrayIndex(objectType, arrayName));
420 }
421 void SetObjectArrayStatus(int objectType, int arrayIndex, int status);
422 void SetObjectArrayStatus(int objectType, const char* arrayName, int status);
424
426
432 int GetNumberOfObjectAttributes(int objectType, int objectIndex);
433 const char* GetObjectAttributeName(int objectType, int objectIndex, int attribIndex);
434 int GetObjectAttributeIndex(int objectType, int objectIndex, const char* attribName);
435 int GetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex);
436 int GetObjectAttributeStatus(int objectType, int objectIndex, const char* attribName)
437 {
438 return this->GetObjectAttributeStatus(
439 objectType, objectIndex, this->GetObjectAttributeIndex(objectType, objectIndex, attribName));
440 }
441 void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status);
442 void SetObjectAttributeStatus(int objectType, int objectIndex, const char* attribName, int status)
443 {
444 this->SetObjectAttributeStatus(objectType, objectIndex,
445 this->GetObjectAttributeIndex(objectType, objectIndex, attribName), status);
446 }
448
453
455
461 const char* GetPartArrayName(int arrayIdx);
462 int GetPartArrayID(const char* name);
463 const char* GetPartBlockInfo(int arrayIdx);
464 void SetPartArrayStatus(int index, int flag);
465 void SetPartArrayStatus(const char*, int flag);
466 int GetPartArrayStatus(int index);
467 int GetPartArrayStatus(const char*);
469
471
478 const char* GetMaterialArrayName(int arrayIdx);
479 int GetMaterialArrayID(const char* name);
480 void SetMaterialArrayStatus(int index, int flag);
481 void SetMaterialArrayStatus(const char*, int flag);
483 int GetMaterialArrayStatus(const char*);
485
487
494 const char* GetAssemblyArrayName(int arrayIdx);
495 int GetAssemblyArrayID(const char* name);
496 void SetAssemblyArrayStatus(int index, int flag);
497 void SetAssemblyArrayStatus(const char*, int flag);
499 int GetAssemblyArrayStatus(const char*);
501
503
513 const char* GetHierarchyArrayName(int arrayIdx);
514 void SetHierarchyArrayStatus(int index, int flag);
515 void SetHierarchyArrayStatus(const char*, int flag);
517 int GetHierarchyArrayStatus(const char*);
519
520 vtkGetMacro(DisplayType, int);
521 virtual void SetDisplayType(int type);
522
526 int IsValidVariable(const char* type, const char* name);
527
531 int GetVariableID(const char* type, const char* name);
532
533 void SetAllArrayStatus(int otype, int status);
534 // Helper functions
535 // static int StringsEqual(const char* s1, char* s2);
536 // static void StringUppercase(const char* str, char* upperstr);
537 // static char *StrDupWithNew(const char *s);
538
539 // time series query functions
540 int GetTimeSeriesData(int ID, const char* vName, const char* vType, vtkFloatArray* result);
541
542 int GetNumberOfEdgeBlockArrays() { return this->GetNumberOfObjects(EDGE_BLOCK); }
543 const char* GetEdgeBlockArrayName(int index) { return this->GetObjectName(EDGE_BLOCK, index); }
544 int GetEdgeBlockArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_BLOCK, name); }
545 void SetEdgeBlockArrayStatus(const char* name, int flag)
546 {
547 this->SetObjectStatus(EDGE_BLOCK, name, flag);
548 }
549
550 int GetNumberOfFaceBlockArrays() { return this->GetNumberOfObjects(FACE_BLOCK); }
551 const char* GetFaceBlockArrayName(int index) { return this->GetObjectName(FACE_BLOCK, index); }
552 int GetFaceBlockArrayStatus(const char* name) { return this->GetObjectStatus(FACE_BLOCK, name); }
553 void SetFaceBlockArrayStatus(const char* name, int flag)
554 {
555 this->SetObjectStatus(FACE_BLOCK, name, flag);
556 }
557
558 int GetNumberOfElementBlockArrays() { return this->GetNumberOfObjects(ELEM_BLOCK); }
559 const char* GetElementBlockArrayName(int index) { return this->GetObjectName(ELEM_BLOCK, index); }
560 int GetElementBlockArrayStatus(const char* name)
561 {
562 return this->GetObjectStatus(ELEM_BLOCK, name);
563 }
564 void SetElementBlockArrayStatus(const char* name, int flag)
565 {
566 this->SetObjectStatus(ELEM_BLOCK, name, flag);
567 }
568
569 int GetNumberOfGlobalResultArrays() { return this->GetNumberOfObjectArrays(GLOBAL); }
570 const char* GetGlobalResultArrayName(int index)
571 {
572 return this->GetObjectArrayName(GLOBAL, index);
573 }
574 int GetGlobalResultArrayStatus(const char* name)
575 {
576 return this->GetObjectArrayStatus(GLOBAL, name);
577 }
578 void SetGlobalResultArrayStatus(const char* name, int flag)
579 {
580 this->SetObjectArrayStatus(GLOBAL, name, flag);
581 }
582
583 int GetNumberOfPointResultArrays() { return this->GetNumberOfObjectArrays(NODAL); }
584 const char* GetPointResultArrayName(int index) { return this->GetObjectArrayName(NODAL, index); }
585 int GetPointResultArrayStatus(const char* name)
586 {
587 return this->GetObjectArrayStatus(NODAL, name);
588 }
589 void SetPointResultArrayStatus(const char* name, int flag)
590 {
591 this->SetObjectArrayStatus(NODAL, name, flag);
592 }
593
594 int GetNumberOfEdgeResultArrays() { return this->GetNumberOfObjectArrays(EDGE_BLOCK); }
595 const char* GetEdgeResultArrayName(int index)
596 {
597 return this->GetObjectArrayName(EDGE_BLOCK, index);
598 }
599 int GetEdgeResultArrayStatus(const char* name)
600 {
601 return this->GetObjectArrayStatus(EDGE_BLOCK, name);
602 }
603 void SetEdgeResultArrayStatus(const char* name, int flag)
604 {
605 this->SetObjectArrayStatus(EDGE_BLOCK, name, flag);
606 }
607
608 int GetNumberOfFaceResultArrays() { return this->GetNumberOfObjectArrays(FACE_BLOCK); }
609 const char* GetFaceResultArrayName(int index)
610 {
611 return this->GetObjectArrayName(FACE_BLOCK, index);
612 }
613 int GetFaceResultArrayStatus(const char* name)
614 {
615 return this->GetObjectArrayStatus(FACE_BLOCK, name);
616 }
617 void SetFaceResultArrayStatus(const char* name, int flag)
618 {
619 this->SetObjectArrayStatus(FACE_BLOCK, name, flag);
620 }
621
622 int GetNumberOfElementResultArrays() { return this->GetNumberOfObjectArrays(ELEM_BLOCK); }
623 const char* GetElementResultArrayName(int index)
624 {
625 return this->GetObjectArrayName(ELEM_BLOCK, index);
626 }
627 int GetElementResultArrayStatus(const char* name)
628 {
629 return this->GetObjectArrayStatus(ELEM_BLOCK, name);
630 }
631 void SetElementResultArrayStatus(const char* name, int flag)
632 {
633 this->SetObjectArrayStatus(ELEM_BLOCK, name, flag);
634 }
635
636 int GetNumberOfNodeMapArrays() { return this->GetNumberOfObjects(NODE_MAP); }
637 const char* GetNodeMapArrayName(int index) { return this->GetObjectName(NODE_MAP, index); }
638 int GetNodeMapArrayStatus(const char* name) { return this->GetObjectStatus(NODE_MAP, name); }
639 void SetNodeMapArrayStatus(const char* name, int flag)
640 {
641 this->SetObjectStatus(NODE_MAP, name, flag);
642 }
643
644 int GetNumberOfEdgeMapArrays() { return this->GetNumberOfObjects(EDGE_MAP); }
645 const char* GetEdgeMapArrayName(int index) { return this->GetObjectName(EDGE_MAP, index); }
646 int GetEdgeMapArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_MAP, name); }
647 void SetEdgeMapArrayStatus(const char* name, int flag)
648 {
649 this->SetObjectStatus(EDGE_MAP, name, flag);
650 }
651
652 int GetNumberOfFaceMapArrays() { return this->GetNumberOfObjects(FACE_MAP); }
653 const char* GetFaceMapArrayName(int index) { return this->GetObjectName(FACE_MAP, index); }
654 int GetFaceMapArrayStatus(const char* name) { return this->GetObjectStatus(FACE_MAP, name); }
655 void SetFaceMapArrayStatus(const char* name, int flag)
656 {
657 this->SetObjectStatus(FACE_MAP, name, flag);
658 }
659
660 int GetNumberOfElementMapArrays() { return this->GetNumberOfObjects(ELEM_MAP); }
661 const char* GetElementMapArrayName(int index) { return this->GetObjectName(ELEM_MAP, index); }
662 int GetElementMapArrayStatus(const char* name) { return this->GetObjectStatus(ELEM_MAP, name); }
663 void SetElementMapArrayStatus(const char* name, int flag)
664 {
665 this->SetObjectStatus(ELEM_MAP, name, flag);
666 }
667
668 int GetNumberOfNodeSetArrays() { return this->GetNumberOfObjects(NODE_SET); }
669 const char* GetNodeSetArrayName(int index) { return this->GetObjectName(NODE_SET, index); }
670 int GetNodeSetArrayStatus(const char* name) { return this->GetObjectStatus(NODE_SET, name); }
671 void SetNodeSetArrayStatus(const char* name, int flag)
672 {
673 this->SetObjectStatus(NODE_SET, name, flag);
674 }
675
676 int GetNumberOfSideSetArrays() { return this->GetNumberOfObjects(SIDE_SET); }
677 const char* GetSideSetArrayName(int index) { return this->GetObjectName(SIDE_SET, index); }
678 int GetSideSetArrayStatus(const char* name) { return this->GetObjectStatus(SIDE_SET, name); }
679 void SetSideSetArrayStatus(const char* name, int flag)
680 {
681 this->SetObjectStatus(SIDE_SET, name, flag);
682 }
683
684 int GetNumberOfEdgeSetArrays() { return this->GetNumberOfObjects(EDGE_SET); }
685 const char* GetEdgeSetArrayName(int index) { return this->GetObjectName(EDGE_SET, index); }
686 int GetEdgeSetArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_SET, name); }
687 void SetEdgeSetArrayStatus(const char* name, int flag)
688 {
689 this->SetObjectStatus(EDGE_SET, name, flag);
690 }
691
692 int GetNumberOfFaceSetArrays() { return this->GetNumberOfObjects(FACE_SET); }
693 const char* GetFaceSetArrayName(int index) { return this->GetObjectName(FACE_SET, index); }
694 int GetFaceSetArrayStatus(const char* name) { return this->GetObjectStatus(FACE_SET, name); }
695 void SetFaceSetArrayStatus(const char* name, int flag)
696 {
697 this->SetObjectStatus(FACE_SET, name, flag);
698 }
699
700 int GetNumberOfElementSetArrays() { return this->GetNumberOfObjects(ELEM_SET); }
701 const char* GetElementSetArrayName(int index) { return this->GetObjectName(ELEM_SET, index); }
702 int GetElementSetArrayStatus(const char* name) { return this->GetObjectStatus(ELEM_SET, name); }
703 void SetElementSetArrayStatus(const char* name, int flag)
704 {
705 this->SetObjectStatus(ELEM_SET, name, flag);
706 }
707
708 int GetNumberOfNodeSetResultArrays() { return this->GetNumberOfObjectArrays(NODE_SET); }
709 const char* GetNodeSetResultArrayName(int index)
710 {
711 return this->GetObjectArrayName(NODE_SET, index);
712 }
713 int GetNodeSetResultArrayStatus(const char* name)
714 {
715 return this->GetObjectArrayStatus(NODE_SET, name);
716 }
717 void SetNodeSetResultArrayStatus(const char* name, int flag)
718 {
719 this->SetObjectArrayStatus(NODE_SET, name, flag);
720 }
721
722 int GetNumberOfSideSetResultArrays() { return this->GetNumberOfObjectArrays(SIDE_SET); }
723 const char* GetSideSetResultArrayName(int index)
724 {
725 return this->GetObjectArrayName(SIDE_SET, index);
726 }
727 int GetSideSetResultArrayStatus(const char* name)
728 {
729 return this->GetObjectArrayStatus(SIDE_SET, name);
730 }
731 void SetSideSetResultArrayStatus(const char* name, int flag)
732 {
733 this->SetObjectArrayStatus(SIDE_SET, name, flag);
734 }
735
736 int GetNumberOfEdgeSetResultArrays() { return this->GetNumberOfObjectArrays(EDGE_SET); }
737 const char* GetEdgeSetResultArrayName(int index)
738 {
739 return this->GetObjectArrayName(EDGE_SET, index);
740 }
741 int GetEdgeSetResultArrayStatus(const char* name)
742 {
743 return this->GetObjectArrayStatus(EDGE_SET, name);
744 }
745 void SetEdgeSetResultArrayStatus(const char* name, int flag)
746 {
747 this->SetObjectArrayStatus(EDGE_SET, name, flag);
748 }
749
750 int GetNumberOfFaceSetResultArrays() { return this->GetNumberOfObjectArrays(FACE_SET); }
751 const char* GetFaceSetResultArrayName(int index)
752 {
753 return this->GetObjectArrayName(FACE_SET, index);
754 }
755 int GetFaceSetResultArrayStatus(const char* name)
756 {
757 return this->GetObjectArrayStatus(FACE_SET, name);
758 }
759 void SetFaceSetResultArrayStatus(const char* name, int flag)
760 {
761 this->SetObjectArrayStatus(FACE_SET, name, flag);
762 }
763
764 int GetNumberOfElementSetResultArrays() { return this->GetNumberOfObjectArrays(ELEM_SET); }
765 const char* GetElementSetResultArrayName(int index)
766 {
767 return this->GetObjectArrayName(ELEM_SET, index);
768 }
769 int GetElementSetResultArrayStatus(const char* name)
770 {
771 return this->GetObjectArrayStatus(ELEM_SET, name);
772 }
773 void SetElementSetResultArrayStatus(const char* name, int flag)
774 {
775 this->SetObjectArrayStatus(ELEM_SET, name, flag);
776 }
777
786 void Reset();
787
797
802
806 void SetCacheSize(double CacheSize);
807
811 double GetCacheSize();
812
814
826 void SetSqueezePoints(bool sp);
829
830 virtual void Dump();
831
837
839
842 vtkGetMacro(SILUpdateStamp, int);
844
846
852
854
866
868
875 vtkSetMacro(UseLegacyBlockNames, bool);
876 vtkGetMacro(UseLegacyBlockNames, bool);
877 vtkBooleanMacro(UseLegacyBlockNames, bool);
879protected:
882
883 // helper for finding IDs
884 static int GetIDHelper(const char* arrayName, vtkDataSet* data, int localID, int searchType);
885 static int GetGlobalID(const char* arrayName, vtkDataSet* data, int localID, int searchType);
886
888 vtkGetObjectMacro(Metadata, vtkExodusIIReaderPrivate);
889
895
896 // Time query function. Called by ExecuteInformation().
897 // Fills the TimestepValues array.
899
904
909 // int RequestDataOverTime( vtkInformation *, vtkInformationVector **, vtkInformationVector *);
910
911 // Parameters for controlling what is read in.
912 char* FileName;
915 int TimeStepRange[2];
918
919 // Information specific for exodus files.
920
921 // 1=display Block names
922 // 2=display Part names
923 // 3=display Material names
925
926 // Metadata containing a description of the currently open file.
928
930
931 friend class vtkPExodusIIReader;
932
933private:
934 vtkExodusIIReader(const vtkExodusIIReader&) = delete;
935 void operator=(const vtkExodusIIReader&) = delete;
936
937 void AddDisplacements(vtkUnstructuredGrid* output);
938 int ModeShapesRange[2];
939
940 bool UseLegacyBlockNames;
941};
942
943VTK_ABI_NAMESPACE_END
944#endif
abstract superclass for arrays of numeric data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
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...
const char * GetObjectTypeName(int)
const char * GetNodeSetArrayName(int index)
int GetEdgeBlockArrayStatus(const char *name)
int GetFaceResultArrayStatus(const char *name)
virtual void SetFileId(int f)
virtual int CanReadFile(const char *fname)
Determine if the file can be read with this reader.
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)
virtual void SetFileName(const char *fname)
Specify file name of the Exodus file.
int GetMaterialArrayStatus(const char *)
By default all materials are loaded.
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)
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.
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)
virtual void SetXMLFileName(const char *fname)
Specify file name of the xml file.
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.
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:342
a simple class to control print indentation
Definition vtkIndent.h:108
Key for integer values in vtkInformation.
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
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:139
record modification and/or execution time
dataset represents arbitrary combinations of all possible cell types
static vtkInformationIntegerKey * GLOBAL_VARIABLE()
Exodus reader outputs global variables and global temporal variables, together with some other variab...
static vtkInformationIntegerKey * GLOBAL_TEMPORAL_VARIABLE()
Exodus reader outputs global variables and global temporal variables, together with some other variab...
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:333
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:288
#define VTK_FILEPATH