VTK  9.4.20250509
vtkHDFWriter.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
36#ifndef vtkHDFWriter_h
37#define vtkHDFWriter_h
38
39#include "vtkIOHDFModule.h" // For export macro
40#include "vtkWriter.h"
41
42#include <map>
43#include <memory>
44
45VTK_ABI_NAMESPACE_BEGIN
46
47class vtkCellArray;
49class vtkDataSet;
50class vtkPoints;
51class vtkPointSet;
52class vtkPolyData;
59
60typedef int64_t hid_t;
61
62class VTKIOHDF_EXPORT vtkHDFWriter : public vtkWriter
63{
64
65private:
66 vtkHDFWriter(const vtkHDFWriter&) = delete;
67 void operator=(const vtkHDFWriter&) = delete;
68
69public:
70 static vtkHDFWriter* New();
71 vtkTypeMacro(vtkHDFWriter, vtkWriter);
72 void PrintSelf(ostream& os, vtkIndent indent) override;
73
75
79 vtkGetObjectMacro(Controller, vtkMultiProcessController);
81
83
89
91
95 vtkSetMacro(Overwrite, bool);
96 vtkGetMacro(Overwrite, bool);
98
100
105 vtkSetMacro(WriteAllTimeSteps, bool);
106 vtkGetMacro(WriteAllTimeSteps, bool);
108
110
122 vtkSetMacro(ChunkSize, int);
123 vtkGetMacro(ChunkSize, int);
125
127
139 vtkSetClampMacro(CompressionLevel, int, 0, 9);
140 vtkGetMacro(CompressionLevel, int);
142
144
151 vtkSetMacro(UseExternalComposite, bool);
152 vtkGetMacro(UseExternalComposite, bool);
154
156
170 vtkSetMacro(UseExternalTimeSteps, bool);
171 vtkGetMacro(UseExternalTimeSteps, bool);
173
175
184 vtkSetMacro(UseExternalPartitions, bool);
185 vtkGetMacro(UseExternalPartitions, bool);
187
188protected:
195 vtkInformationVector* outputVector) override;
196
197 int FillInputPortInformation(int port, vtkInformation* info) override;
198
200 vtkInformationVector* outputVector);
201
202 virtual int RequestUpdateExtent(vtkInformation* request, vtkInformationVector** inputVector,
203 vtkInformationVector* outputVector);
204
206 vtkInformationVector* outputVector) override;
207
209 ~vtkHDFWriter() override;
210
211private:
216 void WriteData() override;
217
222 void DispatchDataObject(hid_t group, vtkDataObject* input, unsigned int partId = 0);
223
228 void WriteDistributedMetafile(vtkDataObject* input);
229
231
235 bool WriteDatasetToFile(hid_t group, vtkPolyData* input, unsigned int partId = 0);
236 bool WriteDatasetToFile(hid_t group, vtkUnstructuredGrid* input, unsigned int partId = 0);
237 bool WriteDatasetToFile(hid_t group, vtkPartitionedDataSet* input);
238 bool WriteDatasetToFile(hid_t group, vtkDataObjectTree* input);
240
242
245 bool UpdateStepsGroup(hid_t group, vtkUnstructuredGrid* input);
246 bool UpdateStepsGroup(hid_t group, vtkPolyData* input);
248
250
254 bool InitializeTemporalPolyData(hid_t group);
255 bool InitializeTemporalUnstructuredGrid(hid_t group);
257
259
263 bool InitializeChunkedDatasets(hid_t group, vtkUnstructuredGrid* input);
264 bool InitializeChunkedDatasets(hid_t group, vtkPolyData* input);
265 bool InitializePointDatasets(hid_t group, vtkPoints* input);
266 bool InitializePrimitiveDataset(hid_t group);
268
273 bool AppendNumberOfPoints(hid_t group, vtkPointSet* input);
274
279 bool AppendPoints(hid_t group, vtkPointSet* input);
280
285 bool AppendNumberOfCells(hid_t group, vtkCellArray* input);
286
291 bool AppendNumberOfConnectivityIds(hid_t group, vtkCellArray* input);
292
297 bool AppendCellTypes(hid_t group, vtkUnstructuredGrid* input);
298
303 bool AppendOffsets(hid_t group, vtkCellArray* input);
304
309 bool AppendConnectivity(hid_t group, vtkCellArray* input);
310
315 bool AppendPrimitiveCells(hid_t baseGroup, vtkPolyData* input);
316
318
322 bool AppendDataArrays(hid_t group, vtkDataObject* input, unsigned int partId = 0);
323 bool AppendDataSetAttributes(hid_t group, vtkDataObject* input, unsigned int partId = 0);
324 bool AppendFieldDataArrays(hid_t group, vtkDataObject* input, unsigned int partId = 0);
326
328
332 bool AppendBlocks(hid_t group, vtkPartitionedDataSetCollection* pdc);
334
336
341 bool AppendExternalBlock(vtkDataObject* block, const std::string& blockName);
343
348 bool AppendAssembly(hid_t group, vtkPartitionedDataSetCollection* pdc);
349
355 bool AppendMultiblock(hid_t group, vtkMultiBlockDataSet* mb, int& leafIndex);
356
361 void AppendIterDataObject(vtkDataObjectTreeIterator* treeIter, const int& leafIndex,
362 const std::string& uniqueSubTreeName);
363
370 void AppendCompositeSubfilesDataObject(const std::string& uniqueSubTreeName);
371
373
376 bool AppendDataArrayOffset(hid_t baseGroup, vtkAbstractArray* array, const std::string& arrayName,
377 const std::string& offsetsGroupName);
378 bool AppendDataArraySizeOffset(hid_t baseGroup, vtkAbstractArray* array,
379 const std::string& arrayName, const std::string& offsetsGroupName);
381
385 bool AppendTimeValues(hid_t group);
386
390 bool HasGeometryChangedFromPreviousStep(vtkDataSet* input);
391
395 void UpdatePreviousStepMeshMTime(vtkDataObject* input);
396
397 class Implementation;
398 std::unique_ptr<Implementation> Impl;
399
400 // Configurable properties
401 char* FileName = nullptr;
402 bool Overwrite = true;
403 bool WriteAllTimeSteps = true;
404 bool UseExternalComposite = false;
405 bool UseExternalTimeSteps = false;
406 bool UseExternalPartitions = false;
407 int ChunkSize = 25000;
408 int CompressionLevel = 0;
409
410 // Temporal-related private variables
411 std::vector<double> timeSteps;
412 bool IsTemporal = false;
413 int CurrentTimeIndex = 0;
414 int NumberOfTimeSteps = 1;
415 vtkMTimeType PreviousStepMeshMTime = 0;
416 std::map<vtkIdType, vtkMTimeType> CompositeMeshMTime;
417
418 // Distributed-related variables
419 vtkMultiProcessController* Controller = nullptr;
420 int NbPieces = 1;
421 int CurrentPiece = 0;
422 bool UsesDummyController = false;
423 std::vector<vtkIdType> PointOffsets;
424 std::vector<vtkIdType> CellOffsets;
425 std::vector<vtkIdType> ConnectivityIdOffsets;
426};
427VTK_ABI_NAMESPACE_END
428#endif
Abstract superclass for all arrays.
object to represent cell connectivity
superclass for composite data iterators
provides implementation for most abstract methods in the superclass vtkCompositeDataSet
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
Write a data object to a VTKHDF file.
vtkTypeBool ProcessRequest(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Override vtkWriter's ProcessRequest method, in order to dispatch the request not only to RequestData ...
virtual int RequestUpdateExtent(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkSetFilePathMacro(FileName)
Get/Set the file name of the vtkHDF file.
static vtkHDFWriter * New()
virtual void SetController(vtkMultiProcessController *)
Set and get the controller.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
vtkGetFilePathMacro(FileName)
Get/Set the file name of the vtkHDF file.
~vtkHDFWriter() override
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
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.
Multiprocessing communication superclass.
Composite dataset that groups datasets as a collection.
composite dataset to encapsulates a dataset consisting of partitions.
concrete class for storing a set of points
Definition vtkPointSet.h:98
represent and manipulate 3D points
Definition vtkPoints.h:139
concrete dataset represents vertices, lines, polygons, and triangle strips
dataset represents arbitrary combinations of all possible cell types
abstract class to write data to file(s)
Definition vtkWriter.h:35
virtual void WriteData()=0
int vtkTypeBool
Definition vtkABI.h:64
int64_t hid_t
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287