VTK  9.6.1
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
32
33#ifndef vtkHDFReader_h
34#define vtkHDFReader_h
35
36#include "vtkDataAssembly.h" // For vtkDataAssembly
38#include "vtkIOHDFModule.h" // For export macro
39#include "vtkSmartPointer.h" // For vtkSmartPointer
40
41#include <array> // For storing the time range
42#include <memory> // For std::unique_ptr
43#include <vector> // For storing list of values
44
45VTK_ABI_NAMESPACE_BEGIN
48class vtkCellData;
49class vtkCommand;
52class vtkDataSet;
55class vtkImageData;
56class vtkInformation;
62class vtkPointData;
63class vtkPolyData;
66
72
73class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
74{
75public:
76 static vtkHDFReader* New();
78 void PrintSelf(ostream& os, vtkIndent indent) override;
79
81
87
89
98
106 virtual int CanReadFile(VTK_FILEPATH const char* name);
107
109
115
117
125
127
133
135
139 const char* GetPointArrayName(int index);
140 const char* GetCellArrayName(int index);
142
144
154 vtkGetMacro(Step, vtkIdType);
155 vtkSetMacro(Step, vtkIdType);
156 vtkGetMacro(TimeValue, double);
157 const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
159
161
170 vtkGetMacro(UseCache, bool);
171 vtkSetMacro(UseCache, bool);
172 vtkBooleanMacro(UseCache, bool);
174
176
194 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks or vtkAppendDataSets instead.")
195 vtkGetMacro(MergeParts, bool);
196 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks vtkAppendDataSets instead.")
197 vtkSetMacro(MergeParts, bool);
198 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks vtkAppendDataSets instead.")
199 virtual void MergePartsOn();
200 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks vtkAppendDataSets instead.")
201 virtual void MergePartsOff();
203
205
211 vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
212 vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
214
216
220 void SetAttributeOriginalIdName(vtkIdType attribute, const std::string& name);
222
227
228protected:
230 ~vtkHDFReader() override;
231
235 int CanReadFileVersion(int major, int minor);
236
238
243 int Read(vtkInformation* outInfo, vtkImageData* data);
250 int ReadRecursively(vtkInformation* outInfo, vtkMultiBlockDataSet* data, const std::string& path);
252
253 VTK_DEPRECATED_IN_9_6_0("This method is deprecated, do not use")
254 int Read(const std::vector<vtkIdType>& numberOfPoints,
255 const std::vector<vtkIdType>& numberOfCells,
256 const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
257 vtkIdType startingPointOffset, vtkIdType startingCellOffset,
258 vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
259
264
269 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
270
272
277 vtkInformationVector* outputVector) override;
279 vtkInformationVector* outputVector) override;
281 vtkInformationVector* outputVector) override;
283
288
293
297 char* FileName;
298
303
309
316
321
323
328 double TimeValue = 0.0;
329 std::array<double, 2> TimeRange;
331
336 // VTK_DEPRECATED_IN_9_5_0( )
337 bool MergeParts = false;
338
340
341 bool UseCache = false;
342 struct DataCache;
343 std::shared_ptr<DataCache> Cache;
344
345private:
346 vtkHDFReader(const vtkHDFReader&) = delete;
347 void operator=(const vtkHDFReader&) = delete;
348
349 class Implementation;
350 Implementation* Impl;
351
356 bool ReadData(vtkInformation* outInfo, vtkDataObject* data);
357
363 int Read(const std::vector<vtkIdType>& numberOfPoints,
364 const std::vector<vtkIdType>& numberOfCells,
365 const std::vector<vtkIdType>& numberOfConnectivityIds,
366 const std::vector<vtkIdType>& numberOfFaces,
367 const std::vector<vtkIdType>& numberOfPolyhedronToFaceIds,
368 const std::vector<vtkIdType>& numberOfFaceConnectivityIds,
369 vtkHDFUtilities::TemporalGeometryOffsets& geoOffsets, int filePiece,
370 vtkUnstructuredGrid* pieceData);
371
377 int Read(const std::vector<vtkIdType>& numberOfTrees, const std::vector<vtkIdType>& numberOfCells,
378 const std::vector<vtkIdType>& numberOfDepths, const std::vector<vtkIdType>& descriptorSizes,
379 const vtkHDFUtilities::TemporalHyperTreeGridOffsets& htgTemporalOffsets, int filePiece,
380 vtkHyperTreeGrid* pieceData);
381
385 void SetHasTemporalData(bool useTemporalData);
386
390 void GenerateAssembly();
391
396 bool RetrieveStepsFromAssembly();
397
402 bool RetrieveDataArraysFromAssembly();
403
411 bool AddOriginalIds(vtkDataSetAttributes* attributes, vtkIdType size, const std::string& name);
412
419 void CleanOriginalIds(vtkPartitionedDataSet* output);
420
421 bool MeshGeometryChangedFromPreviousTimeStep = true;
422
424
425 std::map<vtkIdType, std::string> AttributesOriginalIdName{
426 { vtkDataObject::POINT, "__pointsOriginalIds__" },
427 { vtkDataObject::CELL, "__cellOriginalIds__" }, { vtkDataObject::FIELD, "__fieldOriginalIds__" }
428 };
429
430 bool HasTemporalData = false;
431 std::string CompositeCachePath; // Identifier for the current composite piece
432};
433
434VTK_ABI_NAMESPACE_END
435#endif
Abstract superclass for all arrays.
supports function callbacks
represent and manipulate cell attribute data
Definition vtkCellData.h:32
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:57
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:29
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:58
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
Definition vtkPolyData.h:72
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