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
25#include "vtkIOHDFModule.h" // For export macro
26#include <array> // For storing the time range
27#include <vector> // For storing list of values
28
29VTK_ABI_NAMESPACE_BEGIN
32class vtkCommand;
34class vtkDataSet;
36class vtkImageData;
38class vtkInformation;
40class vtkPolyData;
42
57class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
58{
59public:
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
136protected:
138 ~vtkHDFReader() override;
139
144 constexpr static int GetNumberOfAttributeTypes() { return 3; }
145
149 int CanReadFileVersion(int major, int minor);
150
152
156 int Read(vtkInformation* outInfo, vtkImageData* data);
158 int Read(vtkInformation* outInfo, vtkPolyData* data);
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;
192 vtkInformationVector* outputVector) override;
194
199
200private:
201 vtkHDFReader(const vtkHDFReader&) = delete;
202 void operator=(const vtkHDFReader&) = delete;
203
204protected:
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
247VTK_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
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.
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.
const char * GetCellArrayName(int index)
Get the name of the point or cell array with the given index in the input.
virtual vtkDataArraySelection * GetFieldDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
vtkDataSet * GetOutputAsDataSet(int index)
Get the output as a vtkDataSet pointer.
virtual vtkDataArraySelection * GetCellDataArraySelection()
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
const std::array< double, 2 > & GetTimeRange() const
Getters and setters for transient data.
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.
static vtkHDFReader * New()
static constexpr int GetNumberOfAttributeTypes()
How many attribute types we have.
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.
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.
const char * GetPointArrayName(int index)
Get the name of the point or cell array with the given index in the input.
char * FileName
The input file's name.
virtual vtkDataArraySelection * GetPointDataArraySelection()
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.
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
int Read(vtkInformation *outInfo, vtkPolyData *data)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
vtkDataSet * GetOutputAsDataSet()
Get the output as a vtkDataSet pointer.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
topologically and geometrically regular array of data
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:58
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
int vtkIdType
Definition vtkType.h:315
#define VTK_FILEPATH