VTK  9.6.20260206
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
50
51#ifndef vtkHDFReader_h
52#define vtkHDFReader_h
53
54#include "vtkDataAssembly.h" // For vtkDataAssembly
56#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_7_0 VTK_DEPRECATED_IN_9_6_0
57#include "vtkIOHDFModule.h" // For export macro
58#include "vtkSmartPointer.h" // For vtkSmartPointer
59
60#include <array> // For storing the time range
61#include <memory> // For std::unique_ptr
62#include <vector> // For storing list of values
63
64VTK_ABI_NAMESPACE_BEGIN
67class vtkCellData;
68class vtkCommand;
70class vtkDataSet;
73class vtkImageData;
74class vtkInformation;
80class vtkPointData;
81class vtkPolyData;
84
90
91class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
92{
93public:
94 static vtkHDFReader* New();
96 void PrintSelf(ostream& os, vtkIndent indent) override;
97
99
105
107
116
118
127 virtual vtkTypeBool CanReadFile(VTK_FILEPATH const char* name);
128 static bool CanReadFile(vtkResourceStream* stream);
130
132
138
140
148
150
156
158
162 const char* GetPointArrayName(int index);
163 const char* GetCellArrayName(int index);
165
167
177 vtkGetMacro(Step, vtkIdType);
178 vtkSetMacro(Step, vtkIdType);
179 vtkGetMacro(TimeValue, double);
180 const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
182
184
187 VTK_DEPRECATED_IN_9_7_0("Do not use, cache is on by default and should not be disabled.")
188 virtual bool GetUseCache();
189 VTK_DEPRECATED_IN_9_7_0("Do not use, cache is on by default and should not be disabled.")
190 virtual void SetUseCache(bool use);
191 VTK_DEPRECATED_IN_9_7_0("Do not use, cache is on by default and should not be disabled.")
192 virtual void UseCacheOn();
193 VTK_DEPRECATED_IN_9_7_0("Do not use, cache is on by default and should not be disabled.")
194 virtual void UseCacheOff();
196
198
204 vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
205 vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
207
209
212 VTK_DEPRECATED_IN_9_7_0("Do not use, will be removed.")
214 VTK_DEPRECATED_IN_9_7_0("Do not use, will be removed.")
215 void SetAttributeOriginalIdName(vtkIdType attribute, const std::string& name);
217
222
223protected:
225 ~vtkHDFReader() override;
226
230 int CanReadFileVersion(int major, int minor);
231
233
238 int Read(vtkInformation* outInfo, vtkImageData* data);
245 int ReadRecursively(vtkInformation* outInfo, vtkMultiBlockDataSet* data, const std::string& path);
247
248 VTK_DEPRECATED_IN_9_6_0("This method is deprecated, do not use")
249 int Read(const std::vector<vtkIdType>& numberOfPoints,
250 const std::vector<vtkIdType>& numberOfCells,
251 const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
252 vtkIdType startingPointOffset, vtkIdType startingCellOffset,
253 vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
254
258 VTK_DEPRECATED_IN_9_7_0("This method is deprecated, do not use")
260
265 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
266
268
273 vtkInformationVector* outputVector) override;
275 vtkInformationVector* outputVector) override;
277 vtkInformationVector* outputVector) override;
279
284
289
293 char* FileName;
294
299
305
312
317
319
324 double TimeValue = 0.0;
325 std::array<double, 2> TimeRange;
327
329
330 bool UseCache = true;
331 struct DataCache;
332 std::shared_ptr<DataCache> Cache;
333
334private:
335 vtkHDFReader(const vtkHDFReader&) = delete;
336 void operator=(const vtkHDFReader&) = delete;
337
338 class Implementation;
339 Implementation* Impl;
340
345 bool ReadData(vtkInformation* outInfo, vtkDataObject* data);
346
351 bool ReadAMRData(vtkOverlappingAMR* data, unsigned int maxLevel,
352 vtkDataArraySelection* dataArraySelection[3], bool isTemporalData);
353
359 int Read(const std::vector<vtkIdType>& numberOfPoints,
360 const std::vector<vtkIdType>& numberOfCells,
361 const std::vector<vtkIdType>& numberOfConnectivityIds,
362 const std::vector<vtkIdType>& numberOfFaces,
363 const std::vector<vtkIdType>& numberOfPolyhedronToFaceIds,
364 const std::vector<vtkIdType>& numberOfFaceConnectivityIds,
365 vtkHDFUtilities::TemporalGeometryOffsets& geoOffsets, int filePiece,
366 vtkUnstructuredGrid* pieceData);
367
373 int Read(const std::vector<vtkIdType>& numberOfTrees, const std::vector<vtkIdType>& numberOfCells,
374 const std::vector<vtkIdType>& numberOfDepths, const std::vector<vtkIdType>& descriptorSizes,
375 const vtkHDFUtilities::TemporalHyperTreeGridOffsets& htgTemporalOffsets, int filePiece,
376 vtkHyperTreeGrid* pieceData);
377
381 void SetHasTemporalData(bool useTemporalData);
382
386 void GenerateAssembly();
387
392 bool RetrieveStepsFromAssembly();
393
398 bool RetrieveDataArraysFromAssembly();
399
403 bool ReadFieldArrays(vtkDataObject* data);
404
405 vtkSmartPointer<vtkDataObject> OutputCache;
406
407 std::map<vtkIdType, std::string> AttributesOriginalIdName{
408 { vtkDataObject::POINT, "__pointsOriginalIds__" },
409 { vtkDataObject::CELL, "__cellOriginalIds__" }, { vtkDataObject::FIELD, "__fieldOriginalIds__" }
410 };
411
412 bool HasTemporalData = false;
413 std::string CompositeCachePath; // Identifier for the current composite piece
414};
415
416VTK_ABI_NAMESPACE_END
417#endif
Abstract superclass for all arrays.
supports function callbacks
represent and manipulate cell attribute data
superclass for callback/observer methods
Definition vtkCommand.h:384
Store on/off settings for data arrays, etc.
hierarchical representation to use with vtkPartitionedDataSetCollection
general representation of visualization data
represent and manipulate attribute data in a dataset
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
Implementation for the vtkHDFReader.
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.
static bool CanReadFile(vtkResourceStream *stream)
Return true if, after a quick check of file header, it looks like the provided file or stream can be ...
const char * GetCellArrayName(int index)
Get the name of the point or cell array with the given index in the input.
vtkIdType Step
Temporal data properties.
vtkSmartPointer< vtkResourceStream > Stream
The input stream.
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.
std::shared_ptr< DataCache > Cache
virtual vtkDataArraySelection * GetCellDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
virtual bool GetUseCache()
Boolean property determining whether to use the internal cache or not (default is true).
virtual vtkTypeBool CanReadFile(const char *name)
Return true if, after a quick check of file header, it looks like the provided file or stream can be ...
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) for a specialized data type.
virtual void SetUseCache(bool use)
Boolean property determining whether to use the internal cache or not (default is true).
const std::array< double, 2 > & GetTimeRange() const
Getters and setters for temporal data.
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.
vtkMTimeType GetMTime() override
Overridden to take into account mtime from the internal vtkResourceStream.
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()
virtual void UseCacheOff()
Boolean property determining whether to use the internal cache or not (default is true).
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
Assembly used for PartitionedDataSetCollection.
vtkResourceStream * GetStream()
Specify stream to read from When both Stream and Filename are set, stream is used.
unsigned int MaximumLevelsToReadByDefaultForAMR
void SetAttributeOriginalIdName(vtkIdType attribute, const std::string &name)
Get or Set the Original id name of an attribute (POINT, CELL, FIELD...)
vtkIdType NumberOfSteps
Assembly used for PartitionedDataSetCollection.
double TimeValue
Assembly used for PartitionedDataSetCollection.
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.
int ReadRecursively(vtkInformation *outInfo, vtkMultiBlockDataSet *data, const std::string &path)
Reads the 'data' requested in 'outInfo' (through extents or pieces) for a specialized data type.
int SetupInformation(vtkInformation *outInfo)
Setup the information pass in parameter based on current vtkHDF file loaded.
virtual vtkDataArraySelection * GetPointDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
bool GetHasTemporalData()
Getters and setters for temporal data.
vtkDataArraySelection * DataArraySelection[3]
The array selections.
vtkCallbackCommand * SelectionObserver
The observer to modify this object when the array selections are modified.
std::string GetAttributeOriginalIdName(vtkIdType attribute)
Get or Set the Original id name of an attribute (POINT, CELL, FIELD...)
void SetStream(vtkResourceStream *stream)
Specify stream to read from When both Stream and Filename are set, stream is used.
vtkSmartPointer< vtkDataAssembly > Assembly
Assembly used for PartitionedDataSetCollection.
void PrintPieceInformation(vtkInformation *outInfo)
Print update number of pieces, piece number and ghost levels.
vtkDataSet * GetOutputAsDataSet()
Get the output as a vtkDataSet pointer.
virtual void UseCacheOn()
Boolean property determining whether to use the internal cache or not (default is true).
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
A dataset containing a grid of vtkHyperTree instances arranged as a rectilinear grid.
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Composite dataset that organizes datasets into blocks.
a multi-resolution dataset based on vtkCartesianGrid allowing overlaps
Composite dataset that groups datasets as a collection.
composite dataset to encapsulates a dataset consisting of partitions.
represent and manipulate point attribute data
concrete dataset represents vertices, lines, polygons, and triangle strips
Abstract class used for custom streams.
Hold a reference to a vtkObjectBase instance.
dataset represents arbitrary combinations of all possible cell types
Common utility variables and functions for reader and writer of vtkHDF.
int vtkTypeBool
Definition vtkABI.h:64
#define VTK_DEPRECATED_IN_9_7_0(reason)
#define VTK_DEPRECATED_IN_9_6_0(reason)
int vtkIdType
Definition vtkType.h:354
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:309
#define VTK_FILEPATH