VTK  9.2.6
vtkExodusIIReader.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkExodusIIReader.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
15/*----------------------------------------------------------------------------
16 Copyright (c) Sandia Corporation
17 See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18----------------------------------------------------------------------------*/
19
40#ifndef vtkExodusIIReader_h
41#define vtkExodusIIReader_h
42
43#include "vtkIOExodusModule.h" // For export macro
45
46class vtkDataArray;
47class vtkDataSet;
50class vtkFloatArray;
51class vtkGraph;
53class vtkIntArray;
54class vtkPoints;
56
57class VTKIOEXODUS_EXPORT vtkExodusIIReader : public vtkMultiBlockDataSetAlgorithm
58{
59public:
62 void PrintSelf(ostream& os, vtkIndent indent) override;
63
67 virtual int CanReadFile(VTK_FILEPATH const char* fname);
68
69 // virtual void Modified();
70
75
82
84
87 virtual void SetFileName(VTK_FILEPATH const char* fname);
90
92
95 virtual void SetXMLFileName(VTK_FILEPATH const char* fname);
96 vtkGetFilePathMacro(XMLFileName);
98
100
103 vtkSetMacro(TimeStep, int);
104 vtkGetMacro(TimeStep, int);
106
111 void SetModeShape(int val) { this->SetTimeStep(val - 1); }
112
114
120 vtkGetVector2Macro(ModeShapesRange, int);
122
124
129 vtkGetVector2Macro(TimeStepRange, int);
131
133
146 vtkBooleanMacro(GenerateObjectIdCellArray, vtkTypeBool);
147 static const char* GetObjectIdArrayName() { return "ObjectId"; }
149
152 vtkBooleanMacro(GenerateGlobalElementIdArray, vtkTypeBool);
153
156 vtkBooleanMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
157
160 vtkBooleanMacro(GenerateImplicitElementIdArray, vtkTypeBool);
161
164 vtkBooleanMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
165
168 vtkBooleanMacro(GenerateFileIdArray, vtkTypeBool);
169
170 virtual void SetFileId(int f);
172
174
180 enum
181 {
182 SEARCH_TYPE_ELEMENT = 0,
186 ID_NOT_FOUND = -234121312
187 };
188 // NOTE: GetNumberOfObjectTypes must be updated whenever you add an entry to enum ObjectType {...}
190 {
191 // match Exodus macros from exodusII.h and exodusII_ext.h
192 EDGE_BLOCK = 6,
193 FACE_BLOCK = 8,
194 ELEM_BLOCK = 1,
195 NODE_SET = 2,
196 EDGE_SET = 7,
197 FACE_SET = 9,
198 SIDE_SET = 3,
199 ELEM_SET = 10,
200 NODE_MAP = 5,
201 EDGE_MAP = 11,
202 FACE_MAP = 12,
203 ELEM_MAP = 4,
204 GLOBAL = 13,
205 NODAL = 14,
206 // extended values (not in Exodus headers) for use with SetAllArrayStatus:
207 ASSEMBLY = 60,
208 PART = 61,
209 MATERIAL = 62,
210 HIERARCHY = 63,
211 // extended values (not in Exodus headers) for use in cache keys:
212 QA_RECORDS = 103,
213 INFO_RECORDS = 104,
214 GLOBAL_TEMPORAL = 102,
215 NODAL_TEMPORAL = 101,
216 ELEM_BLOCK_TEMPORAL = 100,
217 GLOBAL_CONN = 99,
218 ELEM_BLOCK_ELEM_CONN = 98,
219 ELEM_BLOCK_FACE_CONN =
220 97,
221 ELEM_BLOCK_EDGE_CONN =
222 96,
223 FACE_BLOCK_CONN = 95,
224 EDGE_BLOCK_CONN = 94,
225 ELEM_SET_CONN = 93,
226 SIDE_SET_CONN = 92,
227 FACE_SET_CONN = 91,
228 EDGE_SET_CONN = 90,
229 NODE_SET_CONN = 89,
230 NODAL_COORDS = 88,
231 OBJECT_ID = 87,
232 IMPLICIT_ELEMENT_ID = 108,
233 IMPLICIT_NODE_ID = 107,
234 GLOBAL_ELEMENT_ID =
235 86,
236 GLOBAL_NODE_ID =
237 85,
238 ELEMENT_ID = 84,
239 NODE_ID = 83,
240 NODAL_SQUEEZEMAP = 82,
241 ELEM_BLOCK_ATTRIB = 81,
242 FACE_BLOCK_ATTRIB = 80,
243 EDGE_BLOCK_ATTRIB = 79,
244 FACE_ID = 105,
245 EDGE_ID = 106,
246 ENTITY_COUNTS = 109
247 };
249
250 static const char* GetGlobalElementIdArrayName() { return "GlobalElementId"; }
251 static const char* GetPedigreeElementIdArrayName() { return "PedigreeElementId"; }
252 static int GetGlobalElementID(vtkDataSet* data, int localID);
253 static int GetGlobalElementID(vtkDataSet* data, int localID, int searchType);
254 static const char* GetImplicitElementIdArrayName() { return "ImplicitElementId"; }
255
256 static const char* GetGlobalFaceIdArrayName() { return "GlobalFaceId"; }
257 static const char* GetPedigreeFaceIdArrayName() { return "PedigreeFaceId"; }
258 static int GetGlobalFaceID(vtkDataSet* data, int localID);
259 static int GetGlobalFaceID(vtkDataSet* data, int localID, int searchType);
260 static const char* GetImplicitFaceIdArrayName() { return "ImplicitFaceId"; }
261
262 static const char* GetGlobalEdgeIdArrayName() { return "GlobalEdgeId"; }
263 static const char* GetPedigreeEdgeIdArrayName() { return "PedigreeEdgeId"; }
264 static int GetGlobalEdgeID(vtkDataSet* data, int localID);
265 static int GetGlobalEdgeID(vtkDataSet* data, int localID, int searchType);
266 static const char* GetImplicitEdgeIdArrayName() { return "ImplicitEdgeId"; }
267
269
275 static const char* GetGlobalNodeIdArrayName() { return "GlobalNodeId"; }
276 static const char* GetPedigreeNodeIdArrayName() { return "PedigreeNodeId"; }
277 static int GetGlobalNodeID(vtkDataSet* data, int localID);
278 static int GetGlobalNodeID(vtkDataSet* data, int localID, int searchType);
279 static const char* GetImplicitNodeIdArrayName() { return "ImplicitNodeId"; }
281
286 static const char* GetSideSetSourceElementIdArrayName() { return "SourceElementId"; }
287
292 static const char* GetSideSetSourceElementSideArrayName() { return "SourceElementSide"; }
294
303 vtkBooleanMacro(ApplyDisplacements, vtkTypeBool);
304 virtual void SetDisplacementMagnitude(float s);
307
309
316 vtkBooleanMacro(HasModeShapes, vtkTypeBool);
318
320
327 virtual void SetModeShapeTime(double phase);
330
332
341 vtkBooleanMacro(AnimateModeShapes, vtkTypeBool);
343
345
351 virtual void SetIgnoreFileTime(bool flag);
353 vtkBooleanMacro(IgnoreFileTime, bool);
355
357
360 const char* GetTitle();
364
369
370 int GetObjectTypeFromName(const char* name);
371 const char* GetObjectTypeName(int);
372
374 int GetNumberOfObjects(int objectType);
375 int GetNumberOfEntriesInObject(int objectType, int objectIndex);
376 int GetObjectId(int objectType, int objectIndex);
377 const char* GetObjectName(int objectType, int objectIndex);
378 using Superclass::GetObjectName;
379 int GetObjectIndex(int objectType, const char* objectName);
380 int GetObjectIndex(int objectType, int id);
381 int GetObjectStatus(int objectType, int objectIndex);
382 int GetObjectStatus(int objectType, const char* objectName)
383 {
384 return this->GetObjectStatus(objectType, this->GetObjectIndex(objectType, objectName));
385 }
386 void SetObjectStatus(int objectType, int objectIndex, int status);
387 void SetObjectStatus(int objectType, const char* objectName, int status);
388
390
396 int GetNumberOfObjectArrays(int objectType);
397 const char* GetObjectArrayName(int objectType, int arrayIndex);
398 int GetObjectArrayIndex(int objectType, const char* arrayName);
399 int GetNumberOfObjectArrayComponents(int objectType, int arrayIndex);
400 int GetObjectArrayStatus(int objectType, int arrayIndex);
401 int GetObjectArrayStatus(int objectType, const char* arrayName)
402 {
403 return this->GetObjectArrayStatus(objectType, this->GetObjectArrayIndex(objectType, arrayName));
404 }
405 void SetObjectArrayStatus(int objectType, int arrayIndex, int status);
406 void SetObjectArrayStatus(int objectType, const char* arrayName, int status);
408
410
416 int GetNumberOfObjectAttributes(int objectType, int objectIndex);
417 const char* GetObjectAttributeName(int objectType, int objectIndex, int attribIndex);
418 int GetObjectAttributeIndex(int objectType, int objectIndex, const char* attribName);
419 int GetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex);
420 int GetObjectAttributeStatus(int objectType, int objectIndex, const char* attribName)
421 {
422 return this->GetObjectAttributeStatus(
423 objectType, objectIndex, this->GetObjectAttributeIndex(objectType, objectIndex, attribName));
424 }
425 void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status);
426 void SetObjectAttributeStatus(int objectType, int objectIndex, const char* attribName, int status)
427 {
428 this->SetObjectAttributeStatus(objectType, objectIndex,
429 this->GetObjectAttributeIndex(objectType, objectIndex, attribName), status);
430 }
432
437
439
445 const char* GetPartArrayName(int arrayIdx);
446 int GetPartArrayID(const char* name);
447 const char* GetPartBlockInfo(int arrayIdx);
448 void SetPartArrayStatus(int index, int flag);
449 void SetPartArrayStatus(const char*, int flag);
450 int GetPartArrayStatus(int index);
451 int GetPartArrayStatus(const char*);
453
455
462 const char* GetMaterialArrayName(int arrayIdx);
463 int GetMaterialArrayID(const char* name);
464 void SetMaterialArrayStatus(int index, int flag);
465 void SetMaterialArrayStatus(const char*, int flag);
467 int GetMaterialArrayStatus(const char*);
469
471
478 const char* GetAssemblyArrayName(int arrayIdx);
479 int GetAssemblyArrayID(const char* name);
480 void SetAssemblyArrayStatus(int index, int flag);
481 void SetAssemblyArrayStatus(const char*, int flag);
483 int GetAssemblyArrayStatus(const char*);
485
487
497 const char* GetHierarchyArrayName(int arrayIdx);
498 void SetHierarchyArrayStatus(int index, int flag);
499 void SetHierarchyArrayStatus(const char*, int flag);
501 int GetHierarchyArrayStatus(const char*);
503
504 vtkGetMacro(DisplayType, int);
505 virtual void SetDisplayType(int type);
506
510 int IsValidVariable(const char* type, const char* name);
511
515 int GetVariableID(const char* type, const char* name);
516
517 void SetAllArrayStatus(int otype, int status);
518 // Helper functions
519 // static int StringsEqual(const char* s1, char* s2);
520 // static void StringUppercase(const char* str, char* upperstr);
521 // static char *StrDupWithNew(const char *s);
522
523 // time series query functions
524 int GetTimeSeriesData(int ID, const char* vName, const char* vType, vtkFloatArray* result);
525
526 int GetNumberOfEdgeBlockArrays() { return this->GetNumberOfObjects(EDGE_BLOCK); }
527 const char* GetEdgeBlockArrayName(int index) { return this->GetObjectName(EDGE_BLOCK, index); }
528 int GetEdgeBlockArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_BLOCK, name); }
529 void SetEdgeBlockArrayStatus(const char* name, int flag)
530 {
531 this->SetObjectStatus(EDGE_BLOCK, name, flag);
532 }
533
534 int GetNumberOfFaceBlockArrays() { return this->GetNumberOfObjects(FACE_BLOCK); }
535 const char* GetFaceBlockArrayName(int index) { return this->GetObjectName(FACE_BLOCK, index); }
536 int GetFaceBlockArrayStatus(const char* name) { return this->GetObjectStatus(FACE_BLOCK, name); }
537 void SetFaceBlockArrayStatus(const char* name, int flag)
538 {
539 this->SetObjectStatus(FACE_BLOCK, name, flag);
540 }
541
542 int GetNumberOfElementBlockArrays() { return this->GetNumberOfObjects(ELEM_BLOCK); }
543 const char* GetElementBlockArrayName(int index) { return this->GetObjectName(ELEM_BLOCK, index); }
544 int GetElementBlockArrayStatus(const char* name)
545 {
546 return this->GetObjectStatus(ELEM_BLOCK, name);
547 }
548 void SetElementBlockArrayStatus(const char* name, int flag)
549 {
550 this->SetObjectStatus(ELEM_BLOCK, name, flag);
551 }
552
553 int GetNumberOfGlobalResultArrays() { return this->GetNumberOfObjectArrays(GLOBAL); }
554 const char* GetGlobalResultArrayName(int index)
555 {
556 return this->GetObjectArrayName(GLOBAL, index);
557 }
558 int GetGlobalResultArrayStatus(const char* name)
559 {
560 return this->GetObjectArrayStatus(GLOBAL, name);
561 }
562 void SetGlobalResultArrayStatus(const char* name, int flag)
563 {
564 this->SetObjectArrayStatus(GLOBAL, name, flag);
565 }
566
567 int GetNumberOfPointResultArrays() { return this->GetNumberOfObjectArrays(NODAL); }
568 const char* GetPointResultArrayName(int index) { return this->GetObjectArrayName(NODAL, index); }
569 int GetPointResultArrayStatus(const char* name)
570 {
571 return this->GetObjectArrayStatus(NODAL, name);
572 }
573 void SetPointResultArrayStatus(const char* name, int flag)
574 {
575 this->SetObjectArrayStatus(NODAL, name, flag);
576 }
577
578 int GetNumberOfEdgeResultArrays() { return this->GetNumberOfObjectArrays(EDGE_BLOCK); }
579 const char* GetEdgeResultArrayName(int index)
580 {
581 return this->GetObjectArrayName(EDGE_BLOCK, index);
582 }
583 int GetEdgeResultArrayStatus(const char* name)
584 {
585 return this->GetObjectArrayStatus(EDGE_BLOCK, name);
586 }
587 void SetEdgeResultArrayStatus(const char* name, int flag)
588 {
589 this->SetObjectArrayStatus(EDGE_BLOCK, name, flag);
590 }
591
592 int GetNumberOfFaceResultArrays() { return this->GetNumberOfObjectArrays(FACE_BLOCK); }
593 const char* GetFaceResultArrayName(int index)
594 {
595 return this->GetObjectArrayName(FACE_BLOCK, index);
596 }
597 int GetFaceResultArrayStatus(const char* name)
598 {
599 return this->GetObjectArrayStatus(FACE_BLOCK, name);
600 }
601 void SetFaceResultArrayStatus(const char* name, int flag)
602 {
603 this->SetObjectArrayStatus(FACE_BLOCK, name, flag);
604 }
605
606 int GetNumberOfElementResultArrays() { return this->GetNumberOfObjectArrays(ELEM_BLOCK); }
607 const char* GetElementResultArrayName(int index)
608 {
609 return this->GetObjectArrayName(ELEM_BLOCK, index);
610 }
611 int GetElementResultArrayStatus(const char* name)
612 {
613 return this->GetObjectArrayStatus(ELEM_BLOCK, name);
614 }
615 void SetElementResultArrayStatus(const char* name, int flag)
616 {
617 this->SetObjectArrayStatus(ELEM_BLOCK, name, flag);
618 }
619
620 int GetNumberOfNodeMapArrays() { return this->GetNumberOfObjects(NODE_MAP); }
621 const char* GetNodeMapArrayName(int index) { return this->GetObjectName(NODE_MAP, index); }
622 int GetNodeMapArrayStatus(const char* name) { return this->GetObjectStatus(NODE_MAP, name); }
623 void SetNodeMapArrayStatus(const char* name, int flag)
624 {
625 this->SetObjectStatus(NODE_MAP, name, flag);
626 }
627
628 int GetNumberOfEdgeMapArrays() { return this->GetNumberOfObjects(EDGE_MAP); }
629 const char* GetEdgeMapArrayName(int index) { return this->GetObjectName(EDGE_MAP, index); }
630 int GetEdgeMapArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_MAP, name); }
631 void SetEdgeMapArrayStatus(const char* name, int flag)
632 {
633 this->SetObjectStatus(EDGE_MAP, name, flag);
634 }
635
636 int GetNumberOfFaceMapArrays() { return this->GetNumberOfObjects(FACE_MAP); }
637 const char* GetFaceMapArrayName(int index) { return this->GetObjectName(FACE_MAP, index); }
638 int GetFaceMapArrayStatus(const char* name) { return this->GetObjectStatus(FACE_MAP, name); }
639 void SetFaceMapArrayStatus(const char* name, int flag)
640 {
641 this->SetObjectStatus(FACE_MAP, name, flag);
642 }
643
644 int GetNumberOfElementMapArrays() { return this->GetNumberOfObjects(ELEM_MAP); }
645 const char* GetElementMapArrayName(int index) { return this->GetObjectName(ELEM_MAP, index); }
646 int GetElementMapArrayStatus(const char* name) { return this->GetObjectStatus(ELEM_MAP, name); }
647 void SetElementMapArrayStatus(const char* name, int flag)
648 {
649 this->SetObjectStatus(ELEM_MAP, name, flag);
650 }
651
652 int GetNumberOfNodeSetArrays() { return this->GetNumberOfObjects(NODE_SET); }
653 const char* GetNodeSetArrayName(int index) { return this->GetObjectName(NODE_SET, index); }
654 int GetNodeSetArrayStatus(const char* name) { return this->GetObjectStatus(NODE_SET, name); }
655 void SetNodeSetArrayStatus(const char* name, int flag)
656 {
657 this->SetObjectStatus(NODE_SET, name, flag);
658 }
659
660 int GetNumberOfSideSetArrays() { return this->GetNumberOfObjects(SIDE_SET); }
661 const char* GetSideSetArrayName(int index) { return this->GetObjectName(SIDE_SET, index); }
662 int GetSideSetArrayStatus(const char* name) { return this->GetObjectStatus(SIDE_SET, name); }
663 void SetSideSetArrayStatus(const char* name, int flag)
664 {
665 this->SetObjectStatus(SIDE_SET, name, flag);
666 }
667
668 int GetNumberOfEdgeSetArrays() { return this->GetNumberOfObjects(EDGE_SET); }
669 const char* GetEdgeSetArrayName(int index) { return this->GetObjectName(EDGE_SET, index); }
670 int GetEdgeSetArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_SET, name); }
671 void SetEdgeSetArrayStatus(const char* name, int flag)
672 {
673 this->SetObjectStatus(EDGE_SET, name, flag);
674 }
675
676 int GetNumberOfFaceSetArrays() { return this->GetNumberOfObjects(FACE_SET); }
677 const char* GetFaceSetArrayName(int index) { return this->GetObjectName(FACE_SET, index); }
678 int GetFaceSetArrayStatus(const char* name) { return this->GetObjectStatus(FACE_SET, name); }
679 void SetFaceSetArrayStatus(const char* name, int flag)
680 {
681 this->SetObjectStatus(FACE_SET, name, flag);
682 }
683
684 int GetNumberOfElementSetArrays() { return this->GetNumberOfObjects(ELEM_SET); }
685 const char* GetElementSetArrayName(int index) { return this->GetObjectName(ELEM_SET, index); }
686 int GetElementSetArrayStatus(const char* name) { return this->GetObjectStatus(ELEM_SET, name); }
687 void SetElementSetArrayStatus(const char* name, int flag)
688 {
689 this->SetObjectStatus(ELEM_SET, name, flag);
690 }
691
692 int GetNumberOfNodeSetResultArrays() { return this->GetNumberOfObjectArrays(NODE_SET); }
693 const char* GetNodeSetResultArrayName(int index)
694 {
695 return this->GetObjectArrayName(NODE_SET, index);
696 }
697 int GetNodeSetResultArrayStatus(const char* name)
698 {
699 return this->GetObjectArrayStatus(NODE_SET, name);
700 }
701 void SetNodeSetResultArrayStatus(const char* name, int flag)
702 {
703 this->SetObjectArrayStatus(NODE_SET, name, flag);
704 }
705
706 int GetNumberOfSideSetResultArrays() { return this->GetNumberOfObjectArrays(SIDE_SET); }
707 const char* GetSideSetResultArrayName(int index)
708 {
709 return this->GetObjectArrayName(SIDE_SET, index);
710 }
711 int GetSideSetResultArrayStatus(const char* name)
712 {
713 return this->GetObjectArrayStatus(SIDE_SET, name);
714 }
715 void SetSideSetResultArrayStatus(const char* name, int flag)
716 {
717 this->SetObjectArrayStatus(SIDE_SET, name, flag);
718 }
719
720 int GetNumberOfEdgeSetResultArrays() { return this->GetNumberOfObjectArrays(EDGE_SET); }
721 const char* GetEdgeSetResultArrayName(int index)
722 {
723 return this->GetObjectArrayName(EDGE_SET, index);
724 }
725 int GetEdgeSetResultArrayStatus(const char* name)
726 {
727 return this->GetObjectArrayStatus(EDGE_SET, name);
728 }
729 void SetEdgeSetResultArrayStatus(const char* name, int flag)
730 {
731 this->SetObjectArrayStatus(EDGE_SET, name, flag);
732 }
733
734 int GetNumberOfFaceSetResultArrays() { return this->GetNumberOfObjectArrays(FACE_SET); }
735 const char* GetFaceSetResultArrayName(int index)
736 {
737 return this->GetObjectArrayName(FACE_SET, index);
738 }
739 int GetFaceSetResultArrayStatus(const char* name)
740 {
741 return this->GetObjectArrayStatus(FACE_SET, name);
742 }
743 void SetFaceSetResultArrayStatus(const char* name, int flag)
744 {
745 this->SetObjectArrayStatus(FACE_SET, name, flag);
746 }
747
748 int GetNumberOfElementSetResultArrays() { return this->GetNumberOfObjectArrays(ELEM_SET); }
749 const char* GetElementSetResultArrayName(int index)
750 {
751 return this->GetObjectArrayName(ELEM_SET, index);
752 }
753 int GetElementSetResultArrayStatus(const char* name)
754 {
755 return this->GetObjectArrayStatus(ELEM_SET, name);
756 }
757 void SetElementSetResultArrayStatus(const char* name, int flag)
758 {
759 this->SetObjectArrayStatus(ELEM_SET, name, flag);
760 }
761
770 void Reset();
771
781
786
790 void SetCacheSize(double CacheSize);
791
795 double GetCacheSize();
796
798
810 void SetSqueezePoints(bool sp);
813
814 virtual void Dump();
815
821
823
826 vtkGetMacro(SILUpdateStamp, int);
828
830
836
838
849
851
858 vtkSetMacro(UseLegacyBlockNames, bool);
859 vtkGetMacro(UseLegacyBlockNames, bool);
860 vtkBooleanMacro(UseLegacyBlockNames, bool);
862protected:
865
866 // helper for finding IDs
867 static int GetIDHelper(const char* arrayName, vtkDataSet* data, int localID, int searchType);
868 static int GetGlobalID(const char* arrayName, vtkDataSet* data, int localID, int searchType);
869
871 vtkGetObjectMacro(Metadata, vtkExodusIIReaderPrivate);
872
878
879 // Time query function. Called by ExecuteInformation().
880 // Fills the TimestepValues array.
882
887
892 // int RequestDataOverTime( vtkInformation *, vtkInformationVector **, vtkInformationVector *);
893
894 // Parameters for controlling what is read in.
895 char* FileName;
898 int TimeStepRange[2];
901
902 // Information specific for exodus files.
903
904 // 1=display Block names
905 // 2=display Part names
906 // 3=display Material names
908
909 // Metadata containing a description of the currently open file.
911
913
914 friend class vtkPExodusIIReader;
915
916private:
917 vtkExodusIIReader(const vtkExodusIIReader&) = delete;
918 void operator=(const vtkExodusIIReader&) = delete;
919
920 void AddDisplacements(vtkUnstructuredGrid* output);
921 int ModeShapesRange[2];
922
923 bool UseLegacyBlockNames;
924};
925
926#endif
abstract superclass for arrays of numeric data
abstract class to specify dataset behavior
Definition vtkDataSet.h:63
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:296
a simple class to control print indentation
Definition vtkIndent.h:40
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:46
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:40
record modification and/or execution time
dataset represents arbitrary combinations of all possible cell types
int vtkTypeBool
Definition vtkABI.h:69
int vtkIdType
Definition vtkType.h:332
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287
#define VTK_FILEPATH