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 
69 #ifndef vtkExodusIIReader_h
70 #define vtkExodusIIReader_h
71 
72 #include "vtkIOExodusModule.h" // For export macro
74 
75 class vtkDataArray;
76 class vtkDataSet;
77 class vtkExodusIICache;
79 class vtkFloatArray;
80 class vtkGraph;
82 class vtkIntArray;
83 class vtkPoints;
85 
86 class VTKIOEXODUS_EXPORT vtkExodusIIReader : public vtkMultiBlockDataSetAlgorithm
87 {
88 public:
91  void PrintSelf(ostream& os, vtkIndent indent) override;
92 
96  virtual int CanReadFile(VTK_FILEPATH const char* fname);
97 
98  // virtual void Modified();
99 
103  vtkMTimeType GetMTime() override;
104 
111 
113 
116  virtual void SetFileName(VTK_FILEPATH const char* fname);
119 
121 
124  virtual void SetXMLFileName(VTK_FILEPATH const char* fname);
125  vtkGetFilePathMacro(XMLFileName);
127 
129 
132  vtkSetMacro(TimeStep, int);
133  vtkGetMacro(TimeStep, int);
135 
140  void SetModeShape(int val) { this->SetTimeStep(val - 1); }
141 
143 
149  vtkGetVector2Macro(ModeShapesRange, int);
151 
153 
158  vtkGetVector2Macro(TimeStepRange, int);
160 
162 
175  vtkBooleanMacro(GenerateObjectIdCellArray, vtkTypeBool);
176  static const char* GetObjectIdArrayName() { return "ObjectId"; }
178 
181  vtkBooleanMacro(GenerateGlobalElementIdArray, vtkTypeBool);
182 
185  vtkBooleanMacro(GenerateGlobalNodeIdArray, vtkTypeBool);
186 
189  vtkBooleanMacro(GenerateImplicitElementIdArray, vtkTypeBool);
190 
193  vtkBooleanMacro(GenerateImplicitNodeIdArray, vtkTypeBool);
194 
197  vtkBooleanMacro(GenerateFileIdArray, vtkTypeBool);
198 
199  virtual void SetFileId(int f);
200  int GetFileId();
201 
203 
209  enum
210  {
211  SEARCH_TYPE_ELEMENT = 0,
215  ID_NOT_FOUND = -234121312
216  };
217  // NOTE: GetNumberOfObjectTypes must be updated whenever you add an entry to enum ObjectType {...}
219  {
220  // match Exodus macros from exodusII.h and exodusII_ext.h
221  EDGE_BLOCK = 6,
222  FACE_BLOCK = 8,
223  ELEM_BLOCK = 1,
224  NODE_SET = 2,
225  EDGE_SET = 7,
226  FACE_SET = 9,
227  SIDE_SET = 3,
228  ELEM_SET = 10,
229  NODE_MAP = 5,
230  EDGE_MAP = 11,
231  FACE_MAP = 12,
232  ELEM_MAP = 4,
233  GLOBAL = 13,
234  NODAL = 14,
235  // extended values (not in Exodus headers) for use with SetAllArrayStatus:
236  ASSEMBLY = 60,
237  PART = 61,
238  MATERIAL = 62,
239  HIERARCHY = 63,
240  // extended values (not in Exodus headers) for use in cache keys:
241  QA_RECORDS = 103,
242  INFO_RECORDS = 104,
243  GLOBAL_TEMPORAL = 102,
244  NODAL_TEMPORAL = 101,
245  ELEM_BLOCK_TEMPORAL = 100,
246  GLOBAL_CONN = 99,
247  ELEM_BLOCK_ELEM_CONN = 98,
248  ELEM_BLOCK_FACE_CONN =
249  97,
250  ELEM_BLOCK_EDGE_CONN =
251  96,
252  FACE_BLOCK_CONN = 95,
253  EDGE_BLOCK_CONN = 94,
254  ELEM_SET_CONN = 93,
255  SIDE_SET_CONN = 92,
256  FACE_SET_CONN = 91,
257  EDGE_SET_CONN = 90,
258  NODE_SET_CONN = 89,
259  NODAL_COORDS = 88,
260  OBJECT_ID = 87,
261  IMPLICIT_ELEMENT_ID = 108,
262  IMPLICIT_NODE_ID = 107,
263  GLOBAL_ELEMENT_ID =
264  86,
265  GLOBAL_NODE_ID =
266  85,
267  ELEMENT_ID = 84,
268  NODE_ID = 83,
269  NODAL_SQUEEZEMAP = 82,
270  ELEM_BLOCK_ATTRIB = 81,
271  FACE_BLOCK_ATTRIB = 80,
272  EDGE_BLOCK_ATTRIB = 79,
273  FACE_ID = 105,
274  EDGE_ID = 106,
275  ENTITY_COUNTS = 109
276  };
278 
279  static const char* GetGlobalElementIdArrayName() { return "GlobalElementId"; }
280  static const char* GetPedigreeElementIdArrayName() { return "PedigreeElementId"; }
281  static int GetGlobalElementID(vtkDataSet* data, int localID);
282  static int GetGlobalElementID(vtkDataSet* data, int localID, int searchType);
283  static const char* GetImplicitElementIdArrayName() { return "ImplicitElementId"; }
284 
285  static const char* GetGlobalFaceIdArrayName() { return "GlobalFaceId"; }
286  static const char* GetPedigreeFaceIdArrayName() { return "PedigreeFaceId"; }
287  static int GetGlobalFaceID(vtkDataSet* data, int localID);
288  static int GetGlobalFaceID(vtkDataSet* data, int localID, int searchType);
289  static const char* GetImplicitFaceIdArrayName() { return "ImplicitFaceId"; }
290 
291  static const char* GetGlobalEdgeIdArrayName() { return "GlobalEdgeId"; }
292  static const char* GetPedigreeEdgeIdArrayName() { return "PedigreeEdgeId"; }
293  static int GetGlobalEdgeID(vtkDataSet* data, int localID);
294  static int GetGlobalEdgeID(vtkDataSet* data, int localID, int searchType);
295  static const char* GetImplicitEdgeIdArrayName() { return "ImplicitEdgeId"; }
296 
298 
304  static const char* GetGlobalNodeIdArrayName() { return "GlobalNodeId"; }
305  static const char* GetPedigreeNodeIdArrayName() { return "PedigreeNodeId"; }
306  static int GetGlobalNodeID(vtkDataSet* data, int localID);
307  static int GetGlobalNodeID(vtkDataSet* data, int localID, int searchType);
308  static const char* GetImplicitNodeIdArrayName() { return "ImplicitNodeId"; }
310 
315  static const char* GetSideSetSourceElementIdArrayName() { return "SourceElementId"; }
316 
321  static const char* GetSideSetSourceElementSideArrayName() { return "SourceElementSide"; }
323 
332  vtkBooleanMacro(ApplyDisplacements, vtkTypeBool);
333  virtual void SetDisplacementMagnitude(float s);
336 
338 
343  virtual void SetHasModeShapes(vtkTypeBool ms);
345  vtkBooleanMacro(HasModeShapes, vtkTypeBool);
347 
349 
356  virtual void SetModeShapeTime(double phase);
359 
361 
368  virtual void SetAnimateModeShapes(vtkTypeBool flag);
370  vtkBooleanMacro(AnimateModeShapes, vtkTypeBool);
372 
374 
380  virtual void SetIgnoreFileTime(bool flag);
382  vtkBooleanMacro(IgnoreFileTime, bool);
384 
386 
389  const char* GetTitle();
393 
398 
399  int GetObjectTypeFromName(const char* name);
400  const char* GetObjectTypeName(int);
401 
403  int GetNumberOfObjects(int objectType);
404  int GetNumberOfEntriesInObject(int objectType, int objectIndex);
405  int GetObjectId(int objectType, int objectIndex);
406  const char* GetObjectName(int objectType, int objectIndex);
408  int GetObjectIndex(int objectType, const char* objectName);
409  int GetObjectIndex(int objectType, int id);
410  int GetObjectStatus(int objectType, int objectIndex);
411  int GetObjectStatus(int objectType, const char* objectName)
412  {
413  return this->GetObjectStatus(objectType, this->GetObjectIndex(objectType, objectName));
414  }
415  void SetObjectStatus(int objectType, int objectIndex, int status);
416  void SetObjectStatus(int objectType, const char* objectName, int status);
417 
419 
425  int GetNumberOfObjectArrays(int objectType);
426  const char* GetObjectArrayName(int objectType, int arrayIndex);
427  int GetObjectArrayIndex(int objectType, const char* arrayName);
428  int GetNumberOfObjectArrayComponents(int objectType, int arrayIndex);
429  int GetObjectArrayStatus(int objectType, int arrayIndex);
430  int GetObjectArrayStatus(int objectType, const char* arrayName)
431  {
432  return this->GetObjectArrayStatus(objectType, this->GetObjectArrayIndex(objectType, arrayName));
433  }
434  void SetObjectArrayStatus(int objectType, int arrayIndex, int status);
435  void SetObjectArrayStatus(int objectType, const char* arrayName, int status);
437 
439 
445  int GetNumberOfObjectAttributes(int objectType, int objectIndex);
446  const char* GetObjectAttributeName(int objectType, int objectIndex, int attribIndex);
447  int GetObjectAttributeIndex(int objectType, int objectIndex, const char* attribName);
448  int GetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex);
449  int GetObjectAttributeStatus(int objectType, int objectIndex, const char* attribName)
450  {
451  return this->GetObjectAttributeStatus(
452  objectType, objectIndex, this->GetObjectAttributeIndex(objectType, objectIndex, attribName));
453  }
454  void SetObjectAttributeStatus(int objectType, int objectIndex, int attribIndex, int status);
455  void SetObjectAttributeStatus(int objectType, int objectIndex, const char* attribName, int status)
456  {
457  this->SetObjectAttributeStatus(objectType, objectIndex,
458  this->GetObjectAttributeIndex(objectType, objectIndex, attribName), status);
459  }
461 
466 
468 
474  const char* GetPartArrayName(int arrayIdx);
475  int GetPartArrayID(const char* name);
476  const char* GetPartBlockInfo(int arrayIdx);
477  void SetPartArrayStatus(int index, int flag);
478  void SetPartArrayStatus(const char*, int flag);
480  int GetPartArrayStatus(const char*);
482 
484 
491  const char* GetMaterialArrayName(int arrayIdx);
492  int GetMaterialArrayID(const char* name);
493  void SetMaterialArrayStatus(int index, int flag);
494  void SetMaterialArrayStatus(const char*, int flag);
496  int GetMaterialArrayStatus(const char*);
498 
500 
507  const char* GetAssemblyArrayName(int arrayIdx);
508  int GetAssemblyArrayID(const char* name);
509  void SetAssemblyArrayStatus(int index, int flag);
510  void SetAssemblyArrayStatus(const char*, int flag);
512  int GetAssemblyArrayStatus(const char*);
514 
516 
526  const char* GetHierarchyArrayName(int arrayIdx);
527  void SetHierarchyArrayStatus(int index, int flag);
528  void SetHierarchyArrayStatus(const char*, int flag);
530  int GetHierarchyArrayStatus(const char*);
532 
533  vtkGetMacro(DisplayType, int);
534  virtual void SetDisplayType(int type);
535 
539  int IsValidVariable(const char* type, const char* name);
540 
544  int GetVariableID(const char* type, const char* name);
545 
546  void SetAllArrayStatus(int otype, int status);
547  // Helper functions
548  // static int StringsEqual(const char* s1, char* s2);
549  // static void StringUppercase(const char* str, char* upperstr);
550  // static char *StrDupWithNew(const char *s);
551 
552  // time series query functions
553  int GetTimeSeriesData(int ID, const char* vName, const char* vType, vtkFloatArray* result);
554 
555  int GetNumberOfEdgeBlockArrays() { return this->GetNumberOfObjects(EDGE_BLOCK); }
556  const char* GetEdgeBlockArrayName(int index) { return this->GetObjectName(EDGE_BLOCK, index); }
557  int GetEdgeBlockArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_BLOCK, name); }
558  void SetEdgeBlockArrayStatus(const char* name, int flag)
559  {
560  this->SetObjectStatus(EDGE_BLOCK, name, flag);
561  }
562 
563  int GetNumberOfFaceBlockArrays() { return this->GetNumberOfObjects(FACE_BLOCK); }
564  const char* GetFaceBlockArrayName(int index) { return this->GetObjectName(FACE_BLOCK, index); }
565  int GetFaceBlockArrayStatus(const char* name) { return this->GetObjectStatus(FACE_BLOCK, name); }
566  void SetFaceBlockArrayStatus(const char* name, int flag)
567  {
568  this->SetObjectStatus(FACE_BLOCK, name, flag);
569  }
570 
571  int GetNumberOfElementBlockArrays() { return this->GetNumberOfObjects(ELEM_BLOCK); }
572  const char* GetElementBlockArrayName(int index) { return this->GetObjectName(ELEM_BLOCK, index); }
574  {
575  return this->GetObjectStatus(ELEM_BLOCK, name);
576  }
577  void SetElementBlockArrayStatus(const char* name, int flag)
578  {
579  this->SetObjectStatus(ELEM_BLOCK, name, flag);
580  }
581 
582  int GetNumberOfGlobalResultArrays() { return this->GetNumberOfObjectArrays(GLOBAL); }
584  {
585  return this->GetObjectArrayName(GLOBAL, index);
586  }
588  {
589  return this->GetObjectArrayStatus(GLOBAL, name);
590  }
591  void SetGlobalResultArrayStatus(const char* name, int flag)
592  {
593  this->SetObjectArrayStatus(GLOBAL, name, flag);
594  }
595 
596  int GetNumberOfPointResultArrays() { return this->GetNumberOfObjectArrays(NODAL); }
597  const char* GetPointResultArrayName(int index) { return this->GetObjectArrayName(NODAL, index); }
599  {
600  return this->GetObjectArrayStatus(NODAL, name);
601  }
602  void SetPointResultArrayStatus(const char* name, int flag)
603  {
604  this->SetObjectArrayStatus(NODAL, name, flag);
605  }
606 
607  int GetNumberOfEdgeResultArrays() { return this->GetNumberOfObjectArrays(EDGE_BLOCK); }
608  const char* GetEdgeResultArrayName(int index)
609  {
610  return this->GetObjectArrayName(EDGE_BLOCK, index);
611  }
613  {
614  return this->GetObjectArrayStatus(EDGE_BLOCK, name);
615  }
616  void SetEdgeResultArrayStatus(const char* name, int flag)
617  {
618  this->SetObjectArrayStatus(EDGE_BLOCK, name, flag);
619  }
620 
621  int GetNumberOfFaceResultArrays() { return this->GetNumberOfObjectArrays(FACE_BLOCK); }
622  const char* GetFaceResultArrayName(int index)
623  {
624  return this->GetObjectArrayName(FACE_BLOCK, index);
625  }
627  {
628  return this->GetObjectArrayStatus(FACE_BLOCK, name);
629  }
630  void SetFaceResultArrayStatus(const char* name, int flag)
631  {
632  this->SetObjectArrayStatus(FACE_BLOCK, name, flag);
633  }
634 
635  int GetNumberOfElementResultArrays() { return this->GetNumberOfObjectArrays(ELEM_BLOCK); }
637  {
638  return this->GetObjectArrayName(ELEM_BLOCK, index);
639  }
641  {
642  return this->GetObjectArrayStatus(ELEM_BLOCK, name);
643  }
644  void SetElementResultArrayStatus(const char* name, int flag)
645  {
646  this->SetObjectArrayStatus(ELEM_BLOCK, name, flag);
647  }
648 
649  int GetNumberOfNodeMapArrays() { return this->GetNumberOfObjects(NODE_MAP); }
650  const char* GetNodeMapArrayName(int index) { return this->GetObjectName(NODE_MAP, index); }
651  int GetNodeMapArrayStatus(const char* name) { return this->GetObjectStatus(NODE_MAP, name); }
652  void SetNodeMapArrayStatus(const char* name, int flag)
653  {
654  this->SetObjectStatus(NODE_MAP, name, flag);
655  }
656 
657  int GetNumberOfEdgeMapArrays() { return this->GetNumberOfObjects(EDGE_MAP); }
658  const char* GetEdgeMapArrayName(int index) { return this->GetObjectName(EDGE_MAP, index); }
659  int GetEdgeMapArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_MAP, name); }
660  void SetEdgeMapArrayStatus(const char* name, int flag)
661  {
662  this->SetObjectStatus(EDGE_MAP, name, flag);
663  }
664 
665  int GetNumberOfFaceMapArrays() { return this->GetNumberOfObjects(FACE_MAP); }
666  const char* GetFaceMapArrayName(int index) { return this->GetObjectName(FACE_MAP, index); }
667  int GetFaceMapArrayStatus(const char* name) { return this->GetObjectStatus(FACE_MAP, name); }
668  void SetFaceMapArrayStatus(const char* name, int flag)
669  {
670  this->SetObjectStatus(FACE_MAP, name, flag);
671  }
672 
673  int GetNumberOfElementMapArrays() { return this->GetNumberOfObjects(ELEM_MAP); }
674  const char* GetElementMapArrayName(int index) { return this->GetObjectName(ELEM_MAP, index); }
675  int GetElementMapArrayStatus(const char* name) { return this->GetObjectStatus(ELEM_MAP, name); }
676  void SetElementMapArrayStatus(const char* name, int flag)
677  {
678  this->SetObjectStatus(ELEM_MAP, name, flag);
679  }
680 
681  int GetNumberOfNodeSetArrays() { return this->GetNumberOfObjects(NODE_SET); }
682  const char* GetNodeSetArrayName(int index) { return this->GetObjectName(NODE_SET, index); }
683  int GetNodeSetArrayStatus(const char* name) { return this->GetObjectStatus(NODE_SET, name); }
684  void SetNodeSetArrayStatus(const char* name, int flag)
685  {
686  this->SetObjectStatus(NODE_SET, name, flag);
687  }
688 
689  int GetNumberOfSideSetArrays() { return this->GetNumberOfObjects(SIDE_SET); }
690  const char* GetSideSetArrayName(int index) { return this->GetObjectName(SIDE_SET, index); }
691  int GetSideSetArrayStatus(const char* name) { return this->GetObjectStatus(SIDE_SET, name); }
692  void SetSideSetArrayStatus(const char* name, int flag)
693  {
694  this->SetObjectStatus(SIDE_SET, name, flag);
695  }
696 
697  int GetNumberOfEdgeSetArrays() { return this->GetNumberOfObjects(EDGE_SET); }
698  const char* GetEdgeSetArrayName(int index) { return this->GetObjectName(EDGE_SET, index); }
699  int GetEdgeSetArrayStatus(const char* name) { return this->GetObjectStatus(EDGE_SET, name); }
700  void SetEdgeSetArrayStatus(const char* name, int flag)
701  {
702  this->SetObjectStatus(EDGE_SET, name, flag);
703  }
704 
705  int GetNumberOfFaceSetArrays() { return this->GetNumberOfObjects(FACE_SET); }
706  const char* GetFaceSetArrayName(int index) { return this->GetObjectName(FACE_SET, index); }
707  int GetFaceSetArrayStatus(const char* name) { return this->GetObjectStatus(FACE_SET, name); }
708  void SetFaceSetArrayStatus(const char* name, int flag)
709  {
710  this->SetObjectStatus(FACE_SET, name, flag);
711  }
712 
713  int GetNumberOfElementSetArrays() { return this->GetNumberOfObjects(ELEM_SET); }
714  const char* GetElementSetArrayName(int index) { return this->GetObjectName(ELEM_SET, index); }
715  int GetElementSetArrayStatus(const char* name) { return this->GetObjectStatus(ELEM_SET, name); }
716  void SetElementSetArrayStatus(const char* name, int flag)
717  {
718  this->SetObjectStatus(ELEM_SET, name, flag);
719  }
720 
721  int GetNumberOfNodeSetResultArrays() { return this->GetNumberOfObjectArrays(NODE_SET); }
723  {
724  return this->GetObjectArrayName(NODE_SET, index);
725  }
727  {
728  return this->GetObjectArrayStatus(NODE_SET, name);
729  }
730  void SetNodeSetResultArrayStatus(const char* name, int flag)
731  {
732  this->SetObjectArrayStatus(NODE_SET, name, flag);
733  }
734 
735  int GetNumberOfSideSetResultArrays() { return this->GetNumberOfObjectArrays(SIDE_SET); }
737  {
738  return this->GetObjectArrayName(SIDE_SET, index);
739  }
741  {
742  return this->GetObjectArrayStatus(SIDE_SET, name);
743  }
744  void SetSideSetResultArrayStatus(const char* name, int flag)
745  {
746  this->SetObjectArrayStatus(SIDE_SET, name, flag);
747  }
748 
749  int GetNumberOfEdgeSetResultArrays() { return this->GetNumberOfObjectArrays(EDGE_SET); }
751  {
752  return this->GetObjectArrayName(EDGE_SET, index);
753  }
755  {
756  return this->GetObjectArrayStatus(EDGE_SET, name);
757  }
758  void SetEdgeSetResultArrayStatus(const char* name, int flag)
759  {
760  this->SetObjectArrayStatus(EDGE_SET, name, flag);
761  }
762 
763  int GetNumberOfFaceSetResultArrays() { return this->GetNumberOfObjectArrays(FACE_SET); }
765  {
766  return this->GetObjectArrayName(FACE_SET, index);
767  }
769  {
770  return this->GetObjectArrayStatus(FACE_SET, name);
771  }
772  void SetFaceSetResultArrayStatus(const char* name, int flag)
773  {
774  this->SetObjectArrayStatus(FACE_SET, name, flag);
775  }
776 
777  int GetNumberOfElementSetResultArrays() { return this->GetNumberOfObjectArrays(ELEM_SET); }
779  {
780  return this->GetObjectArrayName(ELEM_SET, index);
781  }
783  {
784  return this->GetObjectArrayStatus(ELEM_SET, name);
785  }
786  void SetElementSetResultArrayStatus(const char* name, int flag)
787  {
788  this->SetObjectArrayStatus(ELEM_SET, name, flag);
789  }
790 
799  void Reset();
800 
810 
814  void ResetCache();
815 
819  void SetCacheSize(double CacheSize);
820 
824  double GetCacheSize();
825 
827 
839  void SetSqueezePoints(bool sp);
842 
843  virtual void Dump();
844 
850 
852 
855  vtkGetMacro(SILUpdateStamp, int);
857 
859 
865 
867 
878 
880 
887  vtkSetMacro(UseLegacyBlockNames, bool);
888  vtkGetMacro(UseLegacyBlockNames, bool);
889  vtkBooleanMacro(UseLegacyBlockNames, bool);
891 protected:
893  ~vtkExodusIIReader() override;
894 
895  // helper for finding IDs
896  static int GetIDHelper(const char* arrayName, vtkDataSet* data, int localID, int searchType);
897  static int GetGlobalID(const char* arrayName, vtkDataSet* data, int localID, int searchType);
898 
900  vtkGetObjectMacro(Metadata, vtkExodusIIReaderPrivate);
901 
906  bool FindXMLFile();
907 
908  // Time query function. Called by ExecuteInformation().
909  // Fills the TimestepValues array.
911 
916 
921  // int RequestDataOverTime( vtkInformation *, vtkInformationVector **, vtkInformationVector *);
922 
923  // Parameters for controlling what is read in.
924  char* FileName;
925  char* XMLFileName;
926  int TimeStep;
927  int TimeStepRange[2];
930 
931  // Information specific for exodus files.
932 
933  // 1=display Block names
934  // 2=display Part names
935  // 3=display Material names
937 
938  // Metadata containing a description of the currently open file.
940 
942 
943  friend class vtkPExodusIIReader;
944 
945 private:
946  vtkExodusIIReader(const vtkExodusIIReader&) = delete;
947  void operator=(const vtkExodusIIReader&) = delete;
948 
949  void AddDisplacements(vtkUnstructuredGrid* output);
950  int ModeShapesRange[2];
951 
952  bool UseLegacyBlockNames;
953 };
954 
955 #endif
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
abstract class to specify dataset behavior
Definition: vtkDataSet.h:172
This class holds metadata for an Exodus file.
Read exodus 2 files .ex2.
virtual void SetGenerateGlobalNodeIdArray(vtkTypeBool g)
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
see vtkAlgorithm for details
const char * GetFaceSetResultArrayName(int index)
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 * GetEdgeMapArrayName(int index)
const char * GetNodeMapArrayName(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 const char * GetImplicitEdgeIdArrayName()
const char * GetEdgeResultArrayName(int index)
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)
static const char * GetPedigreeElementIdArrayName()
int GetNumberOfObjectArrays(int objectType)
By default arrays are not loaded.
int GetEdgeSetResultArrayStatus(const char *name)
void SetElementMapArrayStatus(const char *name, int flag)
void SetElementSetArrayStatus(const char *name, int flag)
const char * GetEdgeSetResultArrayName(int index)
const char * GetSideSetResultArrayName(int index)
void SetSideSetResultArrayStatus(const char *name, int flag)
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)
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()
void SetFaceBlockArrayStatus(const char *name, int flag)
int GetPointResultArrayStatus(const char *name)
int GetPartArrayStatus(int index)
By default all parts are loaded.
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.
const char * GetElementResultArrayName(int index)
static int GetGlobalID(const char *arrayName, vtkDataSet *data, int localID, int searchType)
static const char * GetPedigreeFaceIdArrayName()
vtkTypeBool GetGenerateGlobalElementIdArray()
virtual void SetGenerateImplicitNodeIdArray(vtkTypeBool g)
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.
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 GetNumberOfEdgeSetResultArrays()
int GetElementResultArrayStatus(const char *name)
int GetNumberOfAssemblyArrays()
By default all assemblies are loaded.
const char * GetObjectName(int objectType, int objectIndex)
void SetAssemblyArrayStatus(const char *, int flag)
By default all assemblies are loaded.
int GetNumberOfSideSetResultArrays()
int GetNumberOfNodeSetResultArrays()
int GetEdgeSetArrayStatus(const char *name)
void SetNodeSetResultArrayStatus(const char *name, int flag)
const char * GetGlobalResultArrayName(int index)
void SetMaterialArrayStatus(int index, int flag)
By default all materials are loaded.
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
static vtkInformationIntegerKey * GLOBAL_VARIABLE()
Exodus reader outputs global variables and global temporal variables, together with some other variab...
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.
const char * GetFaceResultArrayName(int index)
vtkTypeBool GetGenerateObjectIdCellArray()
Extra cell data array that can be generated.
int GetNumberOfMaterialArrays()
By default all materials are loaded.
vtkTypeBool GetGenerateImplicitElementIdArray()
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.
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)
static const char * GetGlobalEdgeIdArrayName()
static vtkExodusIIReader * New()
int GetNumberOfEntriesInObject(int objectType, int objectIndex)
const char * GetElementSetArrayName(int index)
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.
const char * GetElementSetResultArrayName(int index)
const char * GetHierarchyArrayName(int arrayIdx)
By default all hierarchy entries are loaded.
vtkTimeStamp FileNameMTime
int GetAssemblyArrayStatus(const char *)
By default all assemblies are loaded.
void SetModeShape(int val)
Convenience method to set the mode-shape which is same as this->SetTimeStep(val-1);.
int GetAssemblyArrayID(const char *name)
By default all assemblies are loaded.
static const char * GetGlobalElementIdArrayName()
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.
const char * GetSideSetArrayName(int index)
const char * GetTitle()
Access to meta data generated by UpdateInformation.
int GetMaterialArrayID(const char *name)
By default all materials 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)
const char * GetFaceMapArrayName(int index)
void SetObjectStatus(int objectType, int objectIndex, int status)
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...
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 vtkInformationIntegerKey * GLOBAL_TEMPORAL_VARIABLE()
Exodus reader outputs global variables and global temporal variables, together with some other variab...
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.
const char * GetFaceBlockArrayName(int index)
double GetModeShapeTime()
Set/Get the time used to animate mode shapes.
vtkTimeStamp XMLFileNameMTime
static const char * GetPedigreeNodeIdArrayName()
Extra point data array that can be generated.
int GetEdgeResultArrayStatus(const char *name)
const char * GetObjectArrayName(int objectType, int arrayIndex)
By default arrays are not loaded.
const char * GetPartArrayName(int arrayIdx)
By default all parts are loaded.
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.
const char * GetEdgeSetArrayName(int index)
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 * GetMaterialArrayName(int arrayIdx)
By default all materials are loaded.
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.
const char * GetElementBlockArrayName(int index)
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.
vtkTypeBool GetHasModeShapes()
Set/Get whether the Exodus sequence number corresponds to time steps or mode shapes.
void SetEdgeSetArrayStatus(const char *name, int flag)
static const char * GetGlobalNodeIdArrayName()
Extra point data array that can be generated.
virtual vtkIdType GetTotalNumberOfFaces()
int GetObjectId(int objectType, int objectIndex)
int GetNumberOfElementSetResultArrays()
const char * GetNodeSetResultArrayName(int index)
int GetNumberOfEdgesInFile()
void SetElementBlockArrayStatus(const char *name, int flag)
virtual void SetModeShapeTime(double phase)
Set/Get the time used to animate mode shapes.
const char * GetFaceSetArrayName(int index)
void SetFaceSetArrayStatus(const char *name, int flag)
void SetGlobalResultArrayStatus(const char *name, int flag)
vtkTypeBool GetGenerateGlobalNodeIdArray()
vtkGetFilePathMacro(XMLFileName)
Specify file name of the xml file.
virtual void SetApplyDisplacements(vtkTypeBool d)
Geometric locations can include displacements.
virtual vtkIdType GetTotalNumberOfElements()
vtkMTimeType GetMTime() override
Return the object's MTime.
void SetHierarchyArrayStatus(const char *, int flag)
By default all hierarchy entries are loaded.
int GetNumberOfObjects(int objectType)
const char * GetObjectAttributeName(int objectType, int objectIndex, int attribIndex)
By default attributes are not loaded.
static const char * GetGlobalFaceIdArrayName()
static const char * GetObjectIdArrayName()
Extra cell data array that can be generated.
static const char * GetImplicitFaceIdArrayName()
virtual void SetIgnoreFileTime(bool flag)
When on, this option ignores the time values assigned to each time step in the file.
const char * GetNodeSetArrayName(int index)
int GetAssemblyArrayStatus(int index)
By default all assemblies are loaded.
int GetObjectTypeFromName(const char *name)
vtkGraph * GetSIL()
SIL describes organization of/relationships between classifications eg.
const char * GetEdgeBlockArrayName(int index)
int GetElementBlockArrayStatus(const char *name)
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.
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()
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...
const char * GetObjectTypeName(int)
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 * GetElementMapArrayName(int index)
const char * GetPartBlockInfo(int arrayIdx)
By default all parts are loaded.
void SetPointResultArrayStatus(const char *name, int flag)
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)
static const char * GetImplicitElementIdArrayName()
void SetSqueezePoints(bool sp)
Should the reader output only points used by elements in the output mesh, or all the points.
static const char * GetPedigreeEdgeIdArrayName()
virtual void SetGenerateFileIdArray(vtkTypeBool f)
void SetFaceSetResultArrayStatus(const char *name, int flag)
int GetNumberOfFaceSetResultArrays()
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.
ObjectType
Extra cell data array that can be generated.
const char * GetAssemblyArrayName(int arrayIdx)
By default all assemblies are loaded.
int GetNumberOfObjectArrayComponents(int objectType, int arrayIndex)
By default arrays are not loaded.
vtkTypeBool GetApplyDisplacements()
Geometric locations can include displacements.
int GetObjectStatus(int objectType, const char *objectName)
int GetNumberOfElementResultArrays()
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)
int GetNodeSetArrayStatus(const char *name)
int GetFaceSetResultArrayStatus(const char *name)
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 * GetPointResultArrayName(int index)
dynamic, self-adjusting array of float
Base class for graph data types.
Definition: vtkGraph.h:345
a simple class to control print indentation
Definition: vtkIndent.h:119
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:155
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:149
record modification and/or execution time
Definition: vtkTimeStamp.h:55
dataset represents arbitrary combinations of all possible cell types
@ type
Definition: vtkX3D.h:522
@ name
Definition: vtkX3D.h:225
@ index
Definition: vtkX3D.h:252
@ data
Definition: vtkX3D.h:321
int vtkTypeBool
Definition: vtkABI.h:69
int vtkIdType
Definition: vtkType.h:332
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:287
#define VTK_FILEPATH