Loading [MathJax]/extensions/tex2jax.js
VTK  9.4.20250403
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
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;
52class vtkImageData;
54class vtkInformation;
59class vtkPointData;
60class vtkPolyData;
62
80class VTKIOHDF_EXPORT vtkHDFReader : public vtkDataObjectAlgorithm
81{
82public:
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 VTK_DEPRECATED_IN_9_4_0("Please use GetTemporalData method instead.")
149 virtual bool GetHasTransientData();
150 bool GetHasTemporalData();
151 vtkGetMacro(NumberOfSteps, vtkIdType);
152 vtkGetMacro(Step, vtkIdType);
153 vtkSetMacro(Step, vtkIdType);
154 vtkGetMacro(TimeValue, double);
155 const std::array<double, 2>& GetTimeRange() const { return this->TimeRange; }
157
159
168 vtkGetMacro(UseCache, bool);
169 vtkSetMacro(UseCache, bool);
170 vtkBooleanMacro(UseCache, bool);
172
174
192 VTK_DEPRECATED_IN_9_5_0("Use vtkMergeBlocks or vtkAppendDataSets instead.")
193 vtkGetMacro(MergeParts, bool);
195 vtkSetMacro(MergeParts, bool);
197 virtual void MergePartsOn();
199 virtual void MergePartsOff();
201
202 vtkSetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
203 vtkGetMacro(MaximumLevelsToReadByDefaultForAMR, unsigned int);
204
206
209 std::string GetAttributeOriginalIdName(vtkIdType attribute);
210 void SetAttributeOriginalIdName(vtkIdType attribute, const std::string& name);
212
213protected:
215 ~vtkHDFReader() override;
216
220 int CanReadFileVersion(int major, int minor);
221
223
227 int Read(vtkInformation* outInfo, vtkImageData* data);
229 int Read(vtkInformation* outInfo, vtkPolyData* data, vtkPartitionedDataSet* pData);
230 int Read(vtkInformation* outInfo, vtkOverlappingAMR* data);
232 int Read(vtkInformation* outInfo, vtkMultiBlockDataSet* data);
233 int ReadRecursively(vtkInformation* outInfo, vtkMultiBlockDataSet* data, const std::string& path);
235
241 int Read(const std::vector<vtkIdType>& numberOfPoints,
242 const std::vector<vtkIdType>& numberOfCells,
243 const std::vector<vtkIdType>& numberOfConnectivityIds, vtkIdType partOffset,
244 vtkIdType startingPointOffset, vtkIdType startingCellOffset,
245 vtkIdType startingConnectctivityIdOffset, int filePiece, vtkUnstructuredGrid* pieceData);
246
250 int AddFieldArrays(vtkDataObject* data);
251
255 static void SelectionModifiedCallback(
256 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
257
259
263 int RequestDataObject(vtkInformation* request, vtkInformationVector** inputVector,
264 vtkInformationVector* outputVector) override;
265 int RequestInformation(vtkInformation* request, vtkInformationVector** inputVector,
266 vtkInformationVector* outputVector) override;
267 int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
268 vtkInformationVector* outputVector) override;
270
274 void PrintPieceInformation(vtkInformation* outInfo);
275
279 int SetupInformation(vtkInformation* outInfo);
280
284 char* FileName;
285
290 vtkDataArraySelection* DataArraySelection[3];
291
296 vtkCallbackCommand* SelectionObserver;
298
301 int WholeExtent[6];
302 double Origin[3];
303 double Spacing[3];
305
310
312
315 // VTK_DEPRECATED_IN_9_4_0( )
316 bool HasTransientData = false;
317 vtkIdType Step = 0;
318 vtkIdType NumberOfSteps = 1;
319 double TimeValue = 0.0;
320 std::array<double, 2> TimeRange;
322
327 // VTK_DEPRECATED_IN_9_5_0( )
328 bool MergeParts = false;
329
330 unsigned int MaximumLevelsToReadByDefaultForAMR = 0;
331
332 bool UseCache = false;
333 struct DataCache;
334 std::shared_ptr<DataCache> Cache;
335
336private:
337 vtkHDFReader(const vtkHDFReader&) = delete;
338 void operator=(const vtkHDFReader&) = delete;
339
340 class Implementation;
341 Implementation* Impl;
342
348 void SetHasTemporalData(bool useTemporalData);
349
353 void GenerateAssembly();
354
359 bool RetrieveStepsFromAssembly();
360
364 void RetrieveDataArraysFromAssembly();
365
373 bool AddOriginalIds(vtkDataSetAttributes* attributes, vtkIdType size, const std::string& name);
374
381 void CleanOriginalIds(vtkPartitionedDataSet* output);
382
383 bool MeshGeometryChangedFromPreviousTimeStep = true;
384
386
387 std::map<vtkIdType, std::string> AttributesOriginalIdName{
388 { vtkDataObject::POINT, "__pointsOriginalIds__" },
389 { vtkDataObject::CELL, "__cellOriginalIds__" }, { vtkDataObject::FIELD, "__fieldOriginalIds__" }
390 };
391
392 bool HasTemporalData = false;
393};
394
395VTK_ABI_NAMESPACE_END
396#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.
Reads data saved using the VTK HDF format which supports all vtkDataSet types (image data,...
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.
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
#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