VTK  9.5.20251218
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
51
52#ifndef vtkHDFReader_h
53#define vtkHDFReader_h
54
55#include "vtkDataAssembly.h" // For vtkDataAssembly
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;
71class vtkDataSet;
74class vtkImageData;
75class vtkInformation;
81class vtkPointData;
82class vtkPolyData;
85
91
92class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
93{
94public:
95 static vtkHDFReader* New();
97 void PrintSelf(ostream& os, vtkIndent indent) override;
98
100
106
108
117
125 virtual int CanReadFile(VTK_FILEPATH const char* name);
126
128
134
136
144
146
152
154
158 const char* GetPointArrayName(int index);
159 const char* GetCellArrayName(int index);
161
163
173 vtkGetMacro(Step, vtkIdType);
174 vtkSetMacro(Step, vtkIdType);
175 vtkGetMacro(TimeValue, double);
176 const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
178
180
189 vtkGetMacro(UseCache, bool);
190 vtkSetMacro(UseCache, bool);
191 vtkBooleanMacro(UseCache, bool);
193
195
213 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks or vtkAppendDataSets instead.")
214 vtkGetMacro(MergeParts, bool);
215 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks vtkAppendDataSets instead.")
216 vtkSetMacro(MergeParts, bool);
217 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks vtkAppendDataSets instead.")
218 virtual void MergePartsOn();
219 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks vtkAppendDataSets instead.")
220 virtual void MergePartsOff();
222
224
230 vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
231 vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
233
235
239 void SetAttributeOriginalIdName(vtkIdType attribute, const std::string& name);
241
246
247protected:
249 ~vtkHDFReader() override;
250
254 int CanReadFileVersion(int major, int minor);
255
257
262 int Read(vtkInformation* outInfo, vtkImageData* data);
269 int ReadRecursively(vtkInformation* outInfo, vtkMultiBlockDataSet* data, const std::string& path);
271
272 VTK_DEPRECATED_IN_9_6_0("This method is deprecated, do not use")
273 int Read(const std::vector<vtkIdType>& numberOfPoints,
274 const std::vector<vtkIdType>& numberOfCells,
275 const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
276 vtkIdType startingPointOffset, vtkIdType startingCellOffset,
277 vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
278
283
288 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
289
291
296 vtkInformationVector* outputVector) override;
298 vtkInformationVector* outputVector) override;
300 vtkInformationVector* outputVector) override;
302
307
312
316 char* FileName;
317
322
328
335
340
342
347 double TimeValue = 0.0;
348 std::array<double, 2> TimeRange;
350
355 // VTK_DEPRECATED_IN_9_5_0( )
356 bool MergeParts = false;
357
359
360 bool UseCache = false;
361 struct DataCache;
362 std::shared_ptr<DataCache> Cache;
363
364private:
365 vtkHDFReader(const vtkHDFReader&) = delete;
366 void operator=(const vtkHDFReader&) = delete;
367
368 class Implementation;
369 Implementation* Impl;
370
375 bool ReadData(vtkInformation* outInfo, vtkDataObject* data);
376
382 int Read(const std::vector<vtkIdType>& numberOfPoints,
383 const std::vector<vtkIdType>& numberOfCells,
384 const std::vector<vtkIdType>& numberOfConnectivityIds,
385 const std::vector<vtkIdType>& numberOfFaces,
386 const std::vector<vtkIdType>& numberOfPolyhedronToFaceIds,
387 const std::vector<vtkIdType>& numberOfFaceConnectivityIds,
388 vtkHDFUtilities::TemporalGeometryOffsets& geoOffsets, int filePiece,
389 vtkUnstructuredGrid* pieceData);
390
396 int Read(const std::vector<vtkIdType>& numberOfTrees, const std::vector<vtkIdType>& numberOfCells,
397 const std::vector<vtkIdType>& numberOfDepths, const std::vector<vtkIdType>& descriptorSizes,
398 const vtkHDFUtilities::TemporalHyperTreeGridOffsets& htgTemporalOffsets, int filePiece,
399 vtkHyperTreeGrid* pieceData);
400
404 void SetHasTemporalData(bool useTemporalData);
405
409 void GenerateAssembly();
410
415 bool RetrieveStepsFromAssembly();
416
421 bool RetrieveDataArraysFromAssembly();
422
430 bool AddOriginalIds(vtkDataSetAttributes* attributes, vtkIdType size, const std::string& name);
431
438 void CleanOriginalIds(vtkPartitionedDataSet* output);
439
440 bool MeshGeometryChangedFromPreviousTimeStep = true;
441
443
444 std::map<vtkIdType, std::string> AttributesOriginalIdName{
445 { vtkDataObject::POINT, "__pointsOriginalIds__" },
446 { vtkDataObject::CELL, "__cellOriginalIds__" }, { vtkDataObject::FIELD, "__fieldOriginalIds__" }
447 };
448
449 bool HasTemporalData = false;
450 std::string CompositeCachePath; // Identifier for the current composite piece
451};
452
453VTK_ABI_NAMESPACE_END
454#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
vtkDataObjectMeshCache is a class to store and reuse the mesh of a vtkDataSet, while forwarding data ...
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.
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 int CanReadFile(const char *name)
Test whether the file (type) with the given name can be read by this reader.
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) for a specialized data type.
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.
virtual void MergePartsOn()
/!
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 MergePartsOff()
/!
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.
bool MergeParts
/!
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.
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.
Allocate and hold a VTK object.
Definition vtkNew.h:167
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.
#define VTK_DEPRECATED_IN_9_6_0(reason)
#define VTK_DEPRECATED_IN_9_5_0(reason)
int vtkIdType
Definition vtkType.h:368
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:323
#define VTK_FILEPATH