VTK  9.5.2
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
30
31#ifndef vtkHDFReader_h
32#define vtkHDFReader_h
33
34#include "vtkDataAssembly.h" // For vtkDataAssembly
36#include "vtkIOHDFModule.h" // For export macro
37#include "vtkSmartPointer.h" // For vtkSmartPointer
38
39#include <array> // For storing the time range
40#include <memory> // For std::unique_ptr
41#include <vector> // For storing list of values
42
43VTK_ABI_NAMESPACE_BEGIN
46class vtkCellData;
47class vtkCommand;
50class vtkDataSet;
53class vtkImageData;
55class vtkInformation;
60class vtkPointData;
61class vtkPolyData;
63
68
69class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
70{
71public:
72 static vtkHDFReader* New();
74 void PrintSelf(ostream& os, vtkIndent indent) override;
75
77
83
91 virtual int CanReadFile(VTK_FILEPATH const char* name);
92
94
100
102
110
112
118
120
124 const char* GetPointArrayName(int index);
125 const char* GetCellArrayName(int index);
127
129
137 VTK_DEPRECATED_IN_9_4_0("Please use GetTemporalData method instead.")
138 virtual bool GetHasTransientData();
141 vtkGetMacro(Step, vtkIdType);
142 vtkSetMacro(Step, vtkIdType);
143 vtkGetMacro(TimeValue, double);
144 const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
146
148
157 vtkGetMacro(UseCache, bool);
158 vtkSetMacro(UseCache, bool);
159 vtkBooleanMacro(UseCache, bool);
161
163
181 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks or vtkAppendDataSets instead.")
182 vtkGetMacro(MergeParts, bool);
183 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks vtkAppendDataSets instead.")
184 vtkSetMacro(MergeParts, bool);
185 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks vtkAppendDataSets instead.")
186 virtual void MergePartsOn();
187 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks vtkAppendDataSets instead.")
188 virtual void MergePartsOff();
190
192
198 vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
199 vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
201
203
207 void SetAttributeOriginalIdName(vtkIdType attribute, const std::string& name);
209
210protected:
212 ~vtkHDFReader() override;
213
217 int CanReadFileVersion(int major, int minor);
218
220
225 int Read(vtkInformation* outInfo, vtkImageData* data);
232 int ReadRecursively(vtkInformation* outInfo, vtkMultiBlockDataSet* data, const std::string& path);
234
240 int Read(const std::vector<vtkIdType>& numberOfPoints,
241 const std::vector<vtkIdType>& numberOfCells,
242 const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
243 vtkIdType startingPointOffset, vtkIdType startingCellOffset,
244 vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
245
250
255 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
256
258
263 vtkInformationVector* outputVector) override;
265 vtkInformationVector* outputVector) override;
267 vtkInformationVector* outputVector) override;
269
274
279
283 char* FileName;
284
290
297
302
304
307 // VTK_DEPRECATED_IN_9_4_0( )
308 bool HasTransientData = false;
311 double TimeValue = 0.0;
312 std::array<double, 2> TimeRange;
314
319 // VTK_DEPRECATED_IN_9_5_0( )
320 bool MergeParts = false;
321
323
324 bool UseCache = false;
325 struct DataCache;
326 std::shared_ptr<DataCache> Cache;
327
328private:
329 vtkHDFReader(const vtkHDFReader&) = delete;
330 void operator=(const vtkHDFReader&) = delete;
331
332 class Implementation;
333 Implementation* Impl;
334
339 bool ReadData(vtkInformation* outInfo, vtkDataObject* data);
340
346 int Read(const std::vector<vtkIdType>& numberOfTrees, const std::vector<vtkIdType>& numberOfCells,
347 const std::vector<vtkIdType>& numberOfDepths, const std::vector<vtkIdType>& descriptorSizes,
348 const vtkHDFUtilities::TemporalHyperTreeGridOffsets& htgTemporalOffsets, int filePiece,
349 vtkHyperTreeGrid* pieceData);
350
356 void SetHasTemporalData(bool useTemporalData);
357
361 void GenerateAssembly();
362
367 bool RetrieveStepsFromAssembly();
368
373 bool RetrieveDataArraysFromAssembly();
374
382 bool AddOriginalIds(vtkDataSetAttributes* attributes, vtkIdType size, const std::string& name);
383
390 void CleanOriginalIds(vtkPartitionedDataSet* output);
391
392 bool MeshGeometryChangedFromPreviousTimeStep = true;
393
395
396 std::map<vtkIdType, std::string> AttributesOriginalIdName{
397 { vtkDataObject::POINT, "__pointsOriginalIds__" },
398 { vtkDataObject::CELL, "__cellOriginalIds__" }, { vtkDataObject::FIELD, "__fieldOriginalIds__" }
399 };
400
401 bool HasTemporalData = false;
402 std::string CompositeCachePath; // Identifier for the current composite piece
403};
404
405VTK_ABI_NAMESPACE_END
406#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:56
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
Assembly used for PartitionedDataSetCollection.
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.
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()
/!
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.
bool HasTransientData
Temporal data properties.
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.
virtual bool GetHasTransientData()
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...).
virtual int CanReadFile(VTK_FILEPATH const char *name)
Test whether the file (type) with the given name can be read by this reader.
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 vtkUniformGrid 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
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_4_0(reason)
#define VTK_DEPRECATED_IN_9_5_0(reason)
int vtkIdType
Definition vtkType.h:332
#define VTK_FILEPATH