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
50#ifndef vtkHDFReader_h
51#define vtkHDFReader_h
52
53#include "vtkDataAssembly.h" // For vtkDataAssembly
55#include "vtkIOHDFModule.h" // For export macro
56#include "vtkSmartPointer.h" // For vtkSmartPointer
57
58#include <array> // For storing the time range
59#include <memory> // For std::unique_ptr
60#include <vector> // For storing list of values
61
62VTK_ABI_NAMESPACE_BEGIN
65class vtkCellData;
66class vtkCommand;
69class vtkDataSet;
72class vtkImageData;
74class vtkInformation;
79class vtkPointData;
80class vtkPolyData;
82
84{
86}
87
88class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
89{
90public:
91 static vtkHDFReader* New();
93 void PrintSelf(ostream& os, vtkIndent indent) override;
94
96
102
110 virtual int CanReadFile(VTK_FILEPATH const char* name);
111
113
119
121
129
131
137
139
143 const char* GetPointArrayName(int index);
144 const char* GetCellArrayName(int index);
146
148
156 VTK_DEPRECATED_IN_9_4_0("Please use GetTemporalData method instead.")
157 virtual bool GetHasTransientData();
158 bool GetHasTemporalData();
159 vtkGetMacro(NumberOfSteps, vtkIdType);
160 vtkGetMacro(Step, vtkIdType);
161 vtkSetMacro(Step, vtkIdType);
162 vtkGetMacro(TimeValue, double);
163 const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
165
167
176 vtkGetMacro(UseCache, bool);
177 vtkSetMacro(UseCache, bool);
178 vtkBooleanMacro(UseCache, bool);
180
182
200 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks or vtkAppendDataSets instead.")
201 vtkGetMacro(MergeParts, bool);
203 vtkSetMacro(MergeParts, bool);
205 virtual void MergePartsOn();
207 virtual void MergePartsOff();
209
211
217 vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
218 vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
220
222
225 std::string GetAttributeOriginalIdName(vtkIdType attribute);
226 void SetAttributeOriginalIdName(vtkIdType attribute, const std::string& name);
228
229protected:
231 ~vtkHDFReader() override;
232
236 int CanReadFileVersion(int major, int minor);
237
239
244 int Read(vtkInformation* outInfo, vtkImageData* data);
246 int Read(vtkInformation* outInfo, vtkPolyData* data, vtkPartitionedDataSet* pData);
247 int Read(vtkInformation* outInfo, vtkHyperTreeGrid* data, vtkPartitionedDataSet* pData);
248 int Read(vtkInformation* outInfo, vtkOverlappingAMR* data);
250 int Read(vtkInformation* outInfo, vtkMultiBlockDataSet* data);
251 int ReadRecursively(vtkInformation* outInfo, vtkMultiBlockDataSet* data, const std::string& path);
253
259 int Read(const std::vector<vtkIdType>& numberOfPoints,
260 const std::vector<vtkIdType>& numberOfCells,
261 const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
262 vtkIdType startingPointOffset, vtkIdType startingCellOffset,
263 vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
264
268 int AddFieldArrays(vtkDataObject* data);
269
273 static void SelectionModifiedCallback(
274 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
275
277
281 int RequestDataObject(vtkInformation* request, vtkInformationVector** inputVector,
282 vtkInformationVector* outputVector) override;
283 int RequestInformation(vtkInformation* request, vtkInformationVector** inputVector,
284 vtkInformationVector* outputVector) override;
285 int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
286 vtkInformationVector* outputVector) override;
288
292 void PrintPieceInformation(vtkInformation* outInfo);
293
297 int SetupInformation(vtkInformation* outInfo);
298
302 char* FileName;
303
308 vtkDataArraySelection* DataArraySelection[3];
309
314 vtkCallbackCommand* SelectionObserver;
316
321
323
326 // VTK_DEPRECATED_IN_9_4_0( )
327 bool HasTransientData = false;
328 vtkIdType Step = 0;
329 vtkIdType NumberOfSteps = 1;
330 double TimeValue = 0.0;
331 std::array<double, 2> TimeRange;
333
338 // VTK_DEPRECATED_IN_9_5_0( )
339 bool MergeParts = false;
340
341 unsigned int MaximumLevelsToReadByDefaultForAMR = 0;
342
343 bool UseCache = false;
344 struct DataCache;
345 std::shared_ptr<DataCache> Cache;
346
347private:
348 vtkHDFReader(const vtkHDFReader&) = delete;
349 void operator=(const vtkHDFReader&) = delete;
350
351 class Implementation;
352 Implementation* Impl;
353
358 bool ReadData(vtkInformation* outInfo, vtkDataObject* data);
359
365 int Read(const std::vector<vtkIdType>& numberOfTrees, const std::vector<vtkIdType>& numberOfCells,
366 const std::vector<vtkIdType>& numberOfDepths, const std::vector<vtkIdType>& descriptorSizes,
367 const vtkHDFUtilities::TemporalHyperTreeGridOffsets& htgTemporalOffsets, int filePiece,
368 vtkHyperTreeGrid* pieceData);
369
375 void SetHasTemporalData(bool useTemporalData);
376
380 void GenerateAssembly();
381
386 bool RetrieveStepsFromAssembly();
387
392 bool RetrieveDataArraysFromAssembly();
393
401 bool AddOriginalIds(vtkDataSetAttributes* attributes, vtkIdType size, const std::string& name);
402
409 void CleanOriginalIds(vtkPartitionedDataSet* output);
410
411 bool MeshGeometryChangedFromPreviousTimeStep = true;
412
414
415 std::map<vtkIdType, std::string> AttributesOriginalIdName{
416 { vtkDataObject::POINT, "__pointsOriginalIds__" },
417 { vtkDataObject::CELL, "__cellOriginalIds__" }, { vtkDataObject::FIELD, "__fieldOriginalIds__" }
418 };
419
420 bool HasTemporalData = false;
421 std::string CompositeCachePath; // Identifier for the current composite piece
422};
423
424VTK_ABI_NAMESPACE_END
425#endif
Abstract superclass for all arrays.
Appends one or more datasets together into a single output vtkPointSet.
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
Superclass for algorithms that produce only data object as output.
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:165
Implementation for the vtkHDFReader.
Read VTK HDF files.
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.
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.
virtual vtkDataArraySelection * GetCellDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
vtkSetFilePathMacro(FileName)
Get/Set the name of the input file.
static vtkHDFReader * New()
vtkGetFilePathMacro(FileName)
Get/Set the name of the input file.
const char * GetPointArrayName(int index)
Get the name of the point or cell array with the given index in the input.
virtual vtkDataArraySelection * GetPointDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
virtual int CanReadFile(VTK_FILEPATH const char *name)
Test whether the file (type) with the given name can be read by this reader.
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.
merges blocks in a composite dataset to a dataset.
Composite dataset that organizes datasets into blocks.
Allocate and hold a VTK object.
Definition vtkNew.h:167
abstract base class for most VTK objects
Definition vtkObject.h:162
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
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