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
21 #ifndef vtkHDFReader_h
22 #define vtkHDFReader_h
23 
24 #include "vtkDataObjectAlgorithm.h"
25 #include "vtkIOHDFModule.h" // For export macro
26 #include <array> // For storing the time range
27 #include <vector> // For storing list of values
28 
29 VTK_ABI_NAMESPACE_BEGIN
30 class vtkAbstractArray;
31 class vtkCallbackCommand;
32 class vtkCommand;
34 class vtkDataSet;
36 class vtkImageData;
38 class vtkInformation;
39 class vtkOverlappingAMR;
40 class vtkPolyData;
42 
57 class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
58 {
59 public:
60  static vtkHDFReader* New();
62  void PrintSelf(ostream& os, vtkIndent indent) override;
63 
65 
71 
79  virtual int CanReadFile(VTK_FILEPATH const char* name);
80 
82 
88 
90 
98 
100 
106 
108 
112  const char* GetPointArrayName(int index);
113  const char* GetCellArrayName(int index);
115 
117 
125  vtkGetMacro(HasTransientData, bool);
126  vtkGetMacro(NumberOfSteps, vtkIdType);
127  vtkGetMacro(Step, vtkIdType);
128  vtkSetMacro(Step, vtkIdType);
129  vtkGetMacro(TimeValue, double);
130  const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
132 
133  vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
134  vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
135 
136 protected:
138  ~vtkHDFReader() override;
139 
144  constexpr static int GetNumberOfAttributeTypes() { return 3; }
145 
149  int CanReadFileVersion(int major, int minor);
150 
152 
161 
166  int Read(const std::vector<vtkIdType>& numberOfPoints,
167  const std::vector<vtkIdType>& numberOfCells,
168  const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
169  vtkIdType startingPointOffset, vtkIdType startingCellOffset,
170  vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
175 
180  vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
181 
183 
188  vtkInformationVector* outputVector) override;
190  vtkInformationVector* outputVector) override;
191  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
192  vtkInformationVector* outputVector) override;
194 
199 
200 private:
201  vtkHDFReader(const vtkHDFReader&) = delete;
202  void operator=(const vtkHDFReader&) = delete;
203 
204 protected:
208  char* FileName;
209 
214  vtkDataArraySelection* DataArraySelection[3];
215 
222 
225  int WholeExtent[6];
226  double Origin[3];
227  double Spacing[3];
229 
231 
234  bool HasTransientData = false;
235  vtkIdType Step = 0;
236  vtkIdType NumberOfSteps = 1;
237  double TimeValue = 0.0;
238  std::array<double, 2> TimeRange;
240 
241  unsigned int MaximumLevelsToReadByDefaultForAMR = 0;
242 
243  class Implementation;
245 };
246 
247 VTK_ABI_NAMESPACE_END
248 #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:64
represent and manipulate attribute data in a dataset
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
Implementation for the vtkHDFReader.
VTKHDF format reader.
Definition: vtkHDFReader.h:58
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:130
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:238
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:208
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:220
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:243
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:144
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:52
a simple class to control print indentation
Definition: vtkIndent.h:38
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
abstract base class for most VTK objects
Definition: vtkObject.h:61
hierarchical dataset of vtkUniformGrids
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:89
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