VTK  9.3.20240329
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
31 #ifndef vtkHDFReader_h
32 #define vtkHDFReader_h
33 
34 #include "vtkDataAssembly.h" // For vtkDataAssembly
35 #include "vtkDataObjectAlgorithm.h"
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 
43 VTK_ABI_NAMESPACE_BEGIN
44 class vtkAbstractArray;
45 class vtkCallbackCommand;
46 class vtkCellData;
47 class vtkCommand;
50 class vtkDataSet;
52 class vtkImageData;
54 class vtkInformation;
56 class vtkOverlappingAMR;
59 class vtkPointData;
60 class vtkPolyData;
62 
80 class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
81 {
82 public:
83  static vtkHDFReader* New();
85  void PrintSelf(ostream& os, vtkIndent indent) override;
86 
88 
94 
102  virtual int CanReadFile(VTK_FILEPATH const char* name);
103 
105 
111 
113 
121 
123 
129 
131 
135  const char* GetPointArrayName(int index);
136  const char* GetCellArrayName(int index);
138 
140 
148  vtkGetMacro(HasTransientData, bool);
149  vtkGetMacro(NumberOfSteps, vtkIdType);
150  vtkGetMacro(Step, vtkIdType);
151  vtkSetMacro(Step, vtkIdType);
152  vtkGetMacro(TimeValue, double);
153  const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
155 
157 
163  vtkGetMacro(UseCache, bool);
164  vtkSetMacro(UseCache, bool);
165  vtkBooleanMacro(UseCache, bool);
167 
169 
180  vtkGetMacro(MergeParts, bool);
181  vtkSetMacro(MergeParts, bool);
182  vtkBooleanMacro(MergeParts, bool);
184 
185  vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
186  vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
187 
189 
195 
196 protected:
198  ~vtkHDFReader() override;
199 
203  int CanReadFileVersion(int major, int minor);
204 
206 
218 
224  int Read(const std::vector<vtkIdType>& numberOfPoints,
225  const std::vector<vtkIdType>& numberOfCells,
226  const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
227  vtkIdType startingPointOffset, vtkIdType startingCellOffset,
228  vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
229 
234 
239  vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
240 
242 
247  vtkInformationVector* outputVector) override;
249  vtkInformationVector* outputVector) override;
250  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
251  vtkInformationVector* outputVector) override;
253 
258 
263 
264 private:
265  vtkHDFReader(const vtkHDFReader&) = delete;
266  void operator=(const vtkHDFReader&) = delete;
267 
268  bool MeshGeometryChangedFromPreviousTimeStep = true;
269 
271  std::map<vtkIdType, std::string> AttributesOriginalIdName{
272  { vtkDataObject::POINT, "__pointsOriginalIds__" },
273  { vtkDataObject::CELL, "__cellOriginalIds__" }, { vtkDataObject::FIELD, "__fieldOriginalIds__" }
274  };
275 
279  void GenerateAssembly();
280 
285  bool RetrieveStepsFromAssembly();
286 
290  void RetrieveDataArraysFromAssembly();
291 
299  bool AddOriginalIds(vtkDataSetAttributes* attributes, vtkIdType size, const std::string& name);
300 
306  void CleanOriginalIds(vtkDataObject* output);
307 
308 protected:
312  char* FileName;
313 
318  vtkDataArraySelection* DataArraySelection[3];
319 
326 
329  int WholeExtent[6];
330  double Origin[3];
331  double Spacing[3];
333 
338 
340 
343  bool HasTransientData = false;
344  vtkIdType Step = 0;
345  vtkIdType NumberOfSteps = 1;
346  double TimeValue = 0.0;
347  std::array<double, 2> TimeRange;
349 
353  bool MergeParts = true;
354 
355  unsigned int MaximumLevelsToReadByDefaultForAMR = 0;
356 
357  class Implementation;
359 
360  bool UseCache = false;
361  struct DataCache;
362  std::shared_ptr<DataCache> Cache;
363 };
364 
365 VTK_ABI_NAMESPACE_END
366 #endif
Abstract superclass for all arrays.
supports function callbacks
represent and manipulate cell attribute data
Definition: vtkCellData.h:141
superclass for callback/observer methods
Definition: vtkCommand.h:384
Store on/off settings for data arrays, etc.
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:166
Implementation for the vtkHDFReader.
Reads data saved using the VTK HDF format which supports all vtkDataSet types (image data,...
Definition: vtkHDFReader.h:81
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.
std::shared_ptr< DataCache > Cache
Definition: vtkHDFReader.h:361
virtual vtkDataArraySelection * GetFieldDataArraySelection()
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).
~vtkHDFReader() override
int Read(const std::vector< vtkIdType > &numberOfPoints, const std::vector< vtkIdType > &numberOfCells, const std::vector< vtkIdType > &numberOfConnectivityIds, vtkIdType partOffset, vtkIdType startingPointOffset, vtkIdType startingCellOffset, vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid *pieceData)
Read 'pieceData' specified by 'filePiece' where number of points, cells and connectivity ids store th...
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.
int Read(vtkInformation *outInfo, vtkPolyData *data, vtkPartitionedDataSet *pData)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Standard functions to specify the type, information and read the data from the file.
vtkDataSet * GetOutputAsDataSet(int index)
Get the output as a vtkDataSet pointer.
vtkGetFilePathMacro(FileName)
Get/Set the name of the input file.
int Read(vtkInformation *outInfo, vtkUnstructuredGrid *data, vtkPartitionedDataSet *pData)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
static void SelectionModifiedCallback(vtkObject *caller, unsigned long eid, void *clientdata, void *calldata)
Modify this object when an array selection is changed.
const std::array< double, 2 > & GetTimeRange() const
Getters and setters for transient data.
Definition: vtkHDFReader.h:153
int Read(vtkInformation *outInfo, vtkMultiBlockDataSet *data)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
int AddFieldArrays(vtkDataObject *data)
Read the field arrays from the file and add them to the dataset.
std::array< double, 2 > TimeRange
Transient data properties.
Definition: vtkHDFReader.h:347
int Read(vtkInformation *outInfo, vtkPartitionedDataSetCollection *data)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
const char * GetCellArrayName(int index)
Get the name of the point or cell array with the given index in the input.
void SetAttributeOriginalIdName(vtkIdType attribute, const std::string &name)
Get or Set the Original id name of an attribute (POINT, CELL, FIELD...)
virtual vtkDataArraySelection * GetPointDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
char * FileName
The input file's name.
Definition: vtkHDFReader.h:312
int ReadRecursively(vtkInformation *outInfo, vtkMultiBlockDataSet *data, const std::string &path)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
int SetupInformation(vtkInformation *outInfo)
Setup the information pass in parameter based on current vtkHDF file loaded.
virtual vtkDataArraySelection * GetCellDataArraySelection()
Get the data array selection tables used to configure which data arrays are loaded by the reader.
vtkCallbackCommand * SelectionObserver
The observer to modify this object when the array selections are modified.
Definition: vtkHDFReader.h:324
vtkDataSet * GetOutputAsDataSet()
Get the output as a vtkDataSet pointer.
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.
int Read(vtkInformation *outInfo, vtkOverlappingAMR *data)
Reads the 'data' requested in 'outInfo' (through extents or pieces).
vtkSmartPointer< vtkDataAssembly > Assembly
Assembly used for PartitionedDataSetCollection.
Definition: vtkHDFReader.h:337
void PrintPieceInformation(vtkInformation *outInfo)
Print update number of pieces, piece number and ghost levels.
Implementation * Impl
Definition: vtkHDFReader.h:357
static vtkHDFReader * New()
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
const char * GetPointArrayName(int index)
Get the name of the point or cell array with the given index in the input.
topologically and geometrically regular array of data
Definition: vtkImageData.h:156
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.
abstract base class for most VTK objects
Definition: vtkObject.h:162
hierarchical dataset of vtkUniformGrids
Composite dataset that groups datasets as a collection.
composite dataset to encapsulates a dataset consisting of partitions.
represent and manipulate point attribute data
Definition: vtkPointData.h:140
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:181
dataset represents arbitrary combinations of all possible cell types
@ name
Definition: vtkX3D.h:219
@ size
Definition: vtkX3D.h:253
@ index
Definition: vtkX3D.h:246
@ data
Definition: vtkX3D.h:315
@ string
Definition: vtkX3D.h:490
int vtkIdType
Definition: vtkType.h:315
#define VTK_FILEPATH