VTK  9.6.20260325
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
117 enum
118 {
121 };
122
124
136 vtkSetClampMacro(PieceDistribution, int, Block, Interleave);
137 vtkGetMacro(PieceDistribution, int);
139
149 virtual vtkTypeBool CanReadFile(VTK_FILEPATH const char* name);
150 static bool CanReadFile(vtkResourceStream* stream);
152
154
160
162
170
172
178
180
184 const char* GetPointArrayName(int index);
185 const char* GetCellArrayName(int index);
187
189
199 vtkGetMacro(Step, vtkIdType);
200 vtkSetMacro(Step, vtkIdType);
201 vtkGetMacro(TimeValue, double);
202 const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
204
206
209 VTK_DEPRECATED_IN_9_7_0("Do not use, cache is on by default and should not be disabled.")
210 virtual bool GetUseCache();
211 VTK_DEPRECATED_IN_9_7_0("Do not use, cache is on by default and should not be disabled.")
212 virtual void SetUseCache(bool use);
213 VTK_DEPRECATED_IN_9_7_0("Do not use, cache is on by default and should not be disabled.")
214 virtual void UseCacheOn();
215 VTK_DEPRECATED_IN_9_7_0("Do not use, cache is on by default and should not be disabled.")
216 virtual void UseCacheOff();
218
220
226 vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
227 vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
229
231
234 VTK_DEPRECATED_IN_9_7_0("Do not use, will be removed.")
236 VTK_DEPRECATED_IN_9_7_0("Do not use, will be removed.")
237 void SetAttributeOriginalIdName(vtkIdType attribute, const std::string& name);
239
244
245protected:
247 ~vtkHDFReader() override;
248
252 int CanReadFileVersion(int major, int minor);
253
255
260 int Read(vtkInformation* outInfo, vtkImageData* data);
267 int ReadRecursively(vtkInformation* outInfo, vtkMultiBlockDataSet* data, const std::string& path);
269
270 VTK_DEPRECATED_IN_9_6_0("This method is deprecated, do not use")
271 int Read(const std::vector<vtkIdType>& numberOfPoints,
272 const std::vector<vtkIdType>& numberOfCells,
273 const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
274 vtkIdType startingPointOffset, vtkIdType startingCellOffset,
275 vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
276
280 VTK_DEPRECATED_IN_9_7_0("This method is deprecated, do not use")
282
287 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
288
290
295 vtkInformationVector* outputVector) override;
297 vtkInformationVector* outputVector) override;
299 vtkInformationVector* outputVector) override;
301
306
311
315 char* FileName;
316
321
327
334
339
341
346 double TimeValue = 0.0;
347 std::array<double, 2> TimeRange;
349
351
352 bool UseCache = true;
353 struct DataCache;
354 std::shared_ptr<DataCache> Cache;
355
356private:
357 vtkHDFReader(const vtkHDFReader&) = delete;
358 void operator=(const vtkHDFReader&) = delete;
359
360 class Implementation;
361 Implementation* Impl;
362
367 bool ReadData(vtkInformation* outInfo, vtkDataObject* data);
368
373 bool ReadAMRData(vtkOverlappingAMR* data, unsigned int maxLevel,
374 vtkDataArraySelection* dataArraySelection[3], bool isTemporalData);
375
381 int Read(const std::vector<vtkIdType>& numberOfPoints,
382 const std::vector<vtkIdType>& numberOfCells,
383 const std::vector<vtkIdType>& numberOfConnectivityIds,
384 const std::vector<vtkIdType>& numberOfFaces,
385 const std::vector<vtkIdType>& numberOfPolyhedronToFaceIds,
386 const std::vector<vtkIdType>& numberOfFaceConnectivityIds,
387 vtkHDFUtilities::TemporalGeometryOffsets& geoOffsets, int filePiece,
388 vtkUnstructuredGrid* pieceData);
389
395 int Read(const std::vector<vtkIdType>& numberOfTrees, const std::vector<vtkIdType>& numberOfCells,
396 const std::vector<vtkIdType>& numberOfDepths, const std::vector<vtkIdType>& descriptorSizes,
397 const vtkHDFUtilities::TemporalHyperTreeGridOffsets& htgTemporalOffsets, int filePiece,
398 vtkHyperTreeGrid* pieceData);
399
403 void SetHasTemporalData(bool useTemporalData);
404
411 std::vector<int> GetPieceAssignmentForDistribution(
412 int pieceIdx, int numDatasets, int numPieces) const;
413
418 void GenerateAssembly();
419
424 bool RetrieveStepsFromAssembly();
425
430 bool RetrieveDataArraysFromAssembly();
431
435 bool ReadFieldArrays(vtkDataObject* data);
436
437 vtkSmartPointer<vtkDataObject> OutputCache;
438
439 std::map<vtkIdType, std::string> AttributesOriginalIdName{
440 { vtkDataObject::POINT, "__pointsOriginalIds__" },
441 { vtkDataObject::CELL, "__cellOriginalIds__" }, { vtkDataObject::FIELD, "__fieldOriginalIds__" }
442 };
443
444 bool HasTemporalData = false;
445 std::string CompositeCachePath; // Identifier for the current composite piece
446 int PieceDistribution = Interleave;
447};
448
449VTK_ABI_NAMESPACE_END
450#endif
Abstract superclass for all arrays.
supports function callbacks
represent and manipulate cell attribute data
superclass for callback/observer methods
Definition vtkCommand.h:385
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)
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:363
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:318
#define VTK_FILEPATH