VTK  9.3.0
vtkHDFReader.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-License-Identifier: BSD-3-Clause
12 #ifndef vtkHDFReader_h
13 #define vtkHDFReader_h
14 
15 #include "vtkDataObjectAlgorithm.h"
16 #include "vtkIOHDFModule.h" // For export macro
17 #include <array> // For storing the time range
18 #include <vector> // For storing list of values
19 
20 VTK_ABI_NAMESPACE_BEGIN
21 class vtkAbstractArray;
22 class vtkCallbackCommand;
23 class vtkCommand;
25 class vtkDataSet;
27 class vtkImageData;
29 class vtkInformation;
30 class vtkOverlappingAMR;
31 class vtkPolyData;
33 
48 class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
49 {
50 public:
51  static vtkHDFReader* New();
53  void PrintSelf(ostream& os, vtkIndent indent) override;
54 
56 
62 
70  virtual int CanReadFile(VTK_FILEPATH const char* name);
71 
73 
79 
81 
89 
91 
97 
99 
103  const char* GetPointArrayName(int index);
104  const char* GetCellArrayName(int index);
106 
108 
116  vtkGetMacro(HasTransientData, bool);
117  vtkGetMacro(NumberOfSteps, vtkIdType);
118  vtkGetMacro(Step, vtkIdType);
119  vtkSetMacro(Step, vtkIdType);
120  vtkGetMacro(TimeValue, double);
121  const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
123 
124  vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
125  vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
126 
127 protected:
129  ~vtkHDFReader() override;
130 
135  constexpr static int GetNumberOfAttributeTypes() { return 3; }
136 
140  int CanReadFileVersion(int major, int minor);
141 
143 
152 
157  int Read(const std::vector<vtkIdType>& numberOfPoints,
158  const std::vector<vtkIdType>& numberOfCells,
159  const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
160  vtkIdType startingPointOffset, vtkIdType startingCellOffset,
161  vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
166 
171  vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
172 
174 
179  vtkInformationVector* outputVector) override;
181  vtkInformationVector* outputVector) override;
182  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
183  vtkInformationVector* outputVector) override;
185 
190 
191 private:
192  vtkHDFReader(const vtkHDFReader&) = delete;
193  void operator=(const vtkHDFReader&) = delete;
194 
195 protected:
199  char* FileName;
200 
205  vtkDataArraySelection* DataArraySelection[3];
206 
213 
216  int WholeExtent[6];
217  double Origin[3];
218  double Spacing[3];
220 
222 
225  bool HasTransientData = false;
226  vtkIdType Step = 0;
227  vtkIdType NumberOfSteps = 1;
228  double TimeValue = 0.0;
229  std::array<double, 2> TimeRange;
231 
232  unsigned int MaximumLevelsToReadByDefaultForAMR = 0;
233 
234  class Implementation;
236 };
237 
238 VTK_ABI_NAMESPACE_END
239 #endif
Abstract superclass for all arrays.
supports function callbacks
superclass for callback/observer methods
Definition: vtkCommand.h:384
Store on/off settings for data arrays, etc.
Superclass for algorithms that produce only data object as output.
general representation of visualization data
Definition: vtkDataObject.h:55
represent and manipulate attribute data in a dataset
abstract class to specify dataset behavior
Definition: vtkDataSet.h:53
Implementation for the vtkHDFReader.
VTKHDF format reader.
Definition: vtkHDFReader.h:49
int GetNumberOfPointArrays()
Get the number of point or cell arrays available in the input.
int GetNumberOfCellArrays()
Get the number of point or cell arrays available in the input.
virtual vtkDataArraySelection * GetFieldDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
int CanReadFileVersion(int major, int minor)
Test if the reader can read a file with the given version number.
int Read(vtkInformation *outInfo, vtkImageData *data)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
~vtkHDFReader() override
int Read(const std::vector< vtkIdType > &numberOfPoints, const std::vector< vtkIdType > &numberOfCells, const std::vector< vtkIdType > &numberOfConnectivityIds, vtkIdType partOffset, vtkIdType startingPointOffset, vtkIdType startingCellOffset, vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid *pieceData)
Read 'pieceData' specified by 'filePiece' where number of points, cells and connectivity ids store th...
vtkSetFilePathMacro(FileName)
Get/Set the name of the input file.
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Standard functions to specify the type, information and read the data from the file.
int RequestDataObject(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Standard functions to specify the type, information and read the data from the file.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Standard functions to specify the type, information and read the data from the file.
vtkDataSet * GetOutputAsDataSet(int index)
Get the output as a vtkDataSet pointer.
vtkGetFilePathMacro(FileName)
Get/Set the name of the input file.
static void SelectionModifiedCallback(vtkObject *caller, unsigned long eid, void *clientdata, void *calldata)
Modify this object when an array selection is changed.
const std::array< double, 2 > & GetTimeRange() const
Getters and setters for transient data.
Definition: vtkHDFReader.h:121
int AddFieldArrays(vtkDataObject *data)
Read the field arrays from the file and add them to the dataset.
std::array< double, 2 > TimeRange
Transient data properties.
Definition: vtkHDFReader.h:229
const char * GetCellArrayName(int index)
Get the name of the point or cell array with the given index in the input.
virtual vtkDataArraySelection * GetPointDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
char * FileName
The input file's name.
Definition: vtkHDFReader.h:199
virtual vtkDataArraySelection * GetCellDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
vtkCallbackCommand * SelectionObserver
The observer to modify this object when the array selections are modified.
Definition: vtkHDFReader.h:211
vtkDataSet * GetOutputAsDataSet()
Get the output as a vtkDataSet pointer.
virtual int CanReadFile(VTK_FILEPATH const char *name)
Test whether the file (type) with the given name can be read by this reader.
int Read(vtkInformation *outInfo, vtkOverlappingAMR *data)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
void PrintPieceInformation(vtkInformation *outInfo)
Print update number of pieces, piece number and ghost levels.
int Read(vtkInformation *outInfo, vtkUnstructuredGrid *data)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
Implementation * Impl
Definition: vtkHDFReader.h:234
int Read(vtkInformation *outInfo, vtkPolyData *data)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
constexpr static int GetNumberOfAttributeTypes()
How many attribute types we have.
Definition: vtkHDFReader.h:135
static vtkHDFReader * New()
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
const char * GetPointArrayName(int index)
Get the name of the point or cell array with the given index in the input.
topologically and geometrically regular array of data
Definition: vtkImageData.h:43
a simple class to control print indentation
Definition: vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
abstract base class for most VTK objects
Definition: vtkObject.h:52
hierarchical dataset of vtkUniformGrids
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:80
dataset represents arbitrary combinations of all possible cell types
@ name
Definition: vtkX3D.h:219
@ index
Definition: vtkX3D.h:246
@ data
Definition: vtkX3D.h:315
int vtkIdType
Definition: vtkType.h:315
#define VTK_FILEPATH