VTK  9.3.20240418
vtkXdmf3DataSet.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
14 #ifndef vtkXdmf3DataSet_h
15 #define vtkXdmf3DataSet_h
16 
17 #include "vtkIOXdmf3Module.h" // For export macro
18 
19 // clang-format off
20 #include "vtk_xdmf3.h"
21 #include VTKXDMF3_HEADER(core/XdmfSharedPtr.hpp)
22 // clang-format on
23 
24 #include <string> //Needed only for XdmfArray::getName :(
25 
26 class XdmfArray;
27 class XdmfAttribute;
28 class XdmfGrid;
29 class XdmfSet;
30 class XdmfTopologyType;
31 class XdmfRegularGrid;
32 class XdmfRectilinearGrid;
33 class XdmfCurvilinearGrid;
34 class XdmfUnstructuredGrid;
35 class XdmfGraph;
36 class XdmfDomain;
37 
38 VTK_ABI_NAMESPACE_BEGIN
41 class vtkDataArray;
42 class vtkDataObject;
43 class vtkDataSet;
44 class vtkImageData;
45 class vtkRectilinearGrid;
46 class vtkStructuredGrid;
48 class vtkPointSet;
50 class vtkDirectedGraph;
51 
52 class VTKIOXDMF3_EXPORT vtkXdmf3DataSet
53 {
54 public:
55  // Common
56 
60  static vtkDataArray* XdmfToVTKArray(XdmfArray* xArray,
61  std::string attrName, // TODO: needed because XdmfArray::getName() misbehaves
62  unsigned int preferredComponents = 0, vtkXdmf3ArrayKeeper* keeper = nullptr);
63 
67  static bool VTKToXdmfArray(
68  vtkDataArray* vArray, XdmfArray* xArray, unsigned int rank = 0, unsigned int* dims = nullptr);
69 
75  vtkXdmf3ArraySelection* cselection, vtkXdmf3ArraySelection* pselection, XdmfGrid* grid,
76  vtkDataObject* dObject, vtkXdmf3ArrayKeeper* keeper = nullptr);
77 
82  static void VTKToXdmfAttributes(vtkDataObject* dObject, XdmfGrid* grid);
83 
85 
88  static unsigned int GetNumberOfPointsPerCell(int vtk_cell_type, bool& fail);
89  static int GetVTKCellType(shared_ptr<const XdmfTopologyType> topologyType);
90  static int GetXdmfCellType(int vtkType);
92 
94 
97  static void SetTime(XdmfGrid* grid, double hasTime, double time);
98  static void SetTime(XdmfGraph* graph, double hasTime, double time);
100 
101  // vtkXdmf3RegularGrid
102 
106  static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
107  vtkXdmf3ArraySelection* pselection, XdmfRegularGrid* grid, vtkImageData* dataSet,
108  vtkXdmf3ArrayKeeper* keeper = nullptr);
109 
113  static void CopyShape(
114  XdmfRegularGrid* grid, vtkImageData* dataSet, vtkXdmf3ArrayKeeper* keeper = nullptr);
115 
119  static void VTKToXdmf(vtkImageData* dataSet, XdmfDomain* domain, bool hasTime, double time,
120  const char* name = nullptr);
121 
122  // vtkXdmf3RectilinearGrid
126  static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
127  vtkXdmf3ArraySelection* pselection, XdmfRectilinearGrid* grid, vtkRectilinearGrid* dataSet,
128  vtkXdmf3ArrayKeeper* keeper = nullptr);
129 
133  static void CopyShape(
134  XdmfRectilinearGrid* grid, vtkRectilinearGrid* dataSet, vtkXdmf3ArrayKeeper* keeper = nullptr);
135 
139  static void VTKToXdmf(vtkRectilinearGrid* dataSet, XdmfDomain* domain, bool hasTime, double time,
140  const char* name = nullptr);
141 
142  // vtkXdmf3CurvilinearGrid
146  static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
147  vtkXdmf3ArraySelection* pselection, XdmfCurvilinearGrid* grid, vtkStructuredGrid* dataSet,
148  vtkXdmf3ArrayKeeper* keeper = nullptr);
149 
153  static void CopyShape(
154  XdmfCurvilinearGrid* grid, vtkStructuredGrid* dataSet, vtkXdmf3ArrayKeeper* keeper = nullptr);
155 
159  static void VTKToXdmf(vtkStructuredGrid* dataSet, XdmfDomain* domain, bool hasTime, double time,
160  const char* name = nullptr);
161 
162  // vtkXdmf3UnstructuredGrid
166  static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
167  vtkXdmf3ArraySelection* pselection, XdmfUnstructuredGrid* grid, vtkUnstructuredGrid* dataSet,
168  vtkXdmf3ArrayKeeper* keeper = nullptr);
169 
173  static void CopyShape(XdmfUnstructuredGrid* grid, vtkUnstructuredGrid* dataSet,
174  vtkXdmf3ArrayKeeper* keeper = nullptr);
175 
179  static void VTKToXdmf(vtkPointSet* dataSet, XdmfDomain* domain, bool hasTime, double time,
180  const char* name = nullptr);
181 
182  // vtkXdmf3Graph
186  static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
187  vtkXdmf3ArraySelection* pselection, XdmfGraph* grid, vtkMutableDirectedGraph* dataSet,
188  vtkXdmf3ArrayKeeper* keeper = nullptr);
189 
193  static void VTKToXdmf(vtkDirectedGraph* dataSet, XdmfDomain* domain, bool hasTime, double time,
194  const char* name = nullptr);
195 
196  // Side Sets
197 
202  static void XdmfToVTKAttributes(
203  /*
204  vtkXdmf3ArraySelection *fselection,
205  vtkXdmf3ArraySelection *cselection,
206  vtkXdmf3ArraySelection *pselection,
207  */
208  XdmfSet* grid, vtkDataObject* dObject, vtkXdmf3ArrayKeeper* keeper = nullptr);
209 
214  static void XdmfSubsetToVTK(XdmfGrid* grid, unsigned int setnum, vtkDataSet* dataSet,
215  vtkUnstructuredGrid* subSet, vtkXdmf3ArrayKeeper* keeper = nullptr);
216 
222  static int GetVTKFiniteElementCellType(unsigned int element_degree,
223  const std::string& element_family, shared_ptr<const XdmfTopologyType> topologyType);
224 
237  shared_ptr<XdmfAttribute> xmfAttribute, vtkDataArray* array, XdmfGrid* grid,
238  vtkXdmf3ArrayKeeper* keeper = nullptr);
239 };
240 
241 VTK_ABI_NAMESPACE_END
242 #endif
243 // VTK-HeaderTest-Exclude: vtkXdmf3DataSet.h
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:155
general representation of visualization data
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
A directed graph.
topologically and geometrically regular array of data
Definition: vtkImageData.h:156
An editable directed graph.
concrete class for storing a set of points
Definition: vtkPointSet.h:98
a dataset that is topologically regular with variable spacing in the three coordinate directions
topologically regular array of data
dataset represents arbitrary combinations of all possible cell types
LRU cache of XDMF Arrays.
helper to identify requested arrays with
dataset level translation between xdmf3 and vtk
static void CopyShape(XdmfCurvilinearGrid *grid, vtkStructuredGrid *dataSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Helper that does topology for XdmfToVTK.
static int GetVTKFiniteElementCellType(unsigned int element_degree, const std::string &element_family, shared_ptr< const XdmfTopologyType > topologyType)
Converts XDMF topology type, finite element family and degree into an equivalent (or approximative) r...
static void VTKToXdmf(vtkImageData *dataSet, XdmfDomain *domain, bool hasTime, double time, const char *name=nullptr)
Populates the Xdmf Grid with the contents of the VTK data set.
static void XdmfToVTK(vtkXdmf3ArraySelection *fselection, vtkXdmf3ArraySelection *cselection, vtkXdmf3ArraySelection *pselection, XdmfRegularGrid *grid, vtkImageData *dataSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Populates the VTK data set with the contents of the Xdmf grid.
static void VTKToXdmf(vtkRectilinearGrid *dataSet, XdmfDomain *domain, bool hasTime, double time, const char *name=nullptr)
Populates the Xdmf Grid with the contents of the VTK data set.
static void VTKToXdmf(vtkPointSet *dataSet, XdmfDomain *domain, bool hasTime, double time, const char *name=nullptr)
Populates the Xdmf Grid with the contents of the VTK data set.
static void SetTime(XdmfGraph *graph, double hasTime, double time)
Helper used in VTKToXdmf to set the time in a Xdmf grid.
static vtkDataArray * XdmfToVTKArray(XdmfArray *xArray, std::string attrName, unsigned int preferredComponents=0, vtkXdmf3ArrayKeeper *keeper=nullptr)
Returns a VTK array corresponding to the Xdmf array it is given.
static unsigned int GetNumberOfPointsPerCell(int vtk_cell_type, bool &fail)
Helpers for Unstructured Grid translation.
static void CopyShape(XdmfRegularGrid *grid, vtkImageData *dataSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Helper that does topology for XdmfToVTK.
static void VTKToXdmf(vtkDirectedGraph *dataSet, XdmfDomain *domain, bool hasTime, double time, const char *name=nullptr)
Populates the Xdmf Grid with the contents of the VTK data set.
static void XdmfSubsetToVTK(XdmfGrid *grid, unsigned int setnum, vtkDataSet *dataSet, vtkUnstructuredGrid *subSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Extracts numbered subset out of grid (grid corresponds to dataSet), and fills in subSet with it.
static void XdmfToVTK(vtkXdmf3ArraySelection *fselection, vtkXdmf3ArraySelection *cselection, vtkXdmf3ArraySelection *pselection, XdmfGraph *grid, vtkMutableDirectedGraph *dataSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Populates the VTK graph with the contents of the Xdmf grid.
static void SetTime(XdmfGrid *grid, double hasTime, double time)
Helper used in VTKToXdmf to set the time in a Xdmf grid.
static void ParseFiniteElementFunction(vtkDataObject *dObject, shared_ptr< XdmfAttribute > xmfAttribute, vtkDataArray *array, XdmfGrid *grid, vtkXdmf3ArrayKeeper *keeper=nullptr)
Parses finite element function defined in Attribute.
static void XdmfToVTKAttributes(vtkXdmf3ArraySelection *fselection, vtkXdmf3ArraySelection *cselection, vtkXdmf3ArraySelection *pselection, XdmfGrid *grid, vtkDataObject *dObject, vtkXdmf3ArrayKeeper *keeper=nullptr)
Populates the given VTK DataObject's attribute arrays with the selected arrays from the Xdmf Grid.
static void XdmfToVTK(vtkXdmf3ArraySelection *fselection, vtkXdmf3ArraySelection *cselection, vtkXdmf3ArraySelection *pselection, XdmfCurvilinearGrid *grid, vtkStructuredGrid *dataSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Populates the VTK data set with the contents of the Xdmf grid.
static int GetVTKCellType(shared_ptr< const XdmfTopologyType > topologyType)
Helpers for Unstructured Grid translation.
static int GetXdmfCellType(int vtkType)
Helpers for Unstructured Grid translation.
static void VTKToXdmf(vtkStructuredGrid *dataSet, XdmfDomain *domain, bool hasTime, double time, const char *name=nullptr)
Populates the Xdmf Grid with the contents of the VTK data set.
static void VTKToXdmfAttributes(vtkDataObject *dObject, XdmfGrid *grid)
Populates the given Xdmf Grid's attribute arrays with the selected arrays from the VTK DataObject.
static void XdmfToVTK(vtkXdmf3ArraySelection *fselection, vtkXdmf3ArraySelection *cselection, vtkXdmf3ArraySelection *pselection, XdmfRectilinearGrid *grid, vtkRectilinearGrid *dataSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Populates the VTK data set with the contents of the Xdmf grid.
static void XdmfToVTK(vtkXdmf3ArraySelection *fselection, vtkXdmf3ArraySelection *cselection, vtkXdmf3ArraySelection *pselection, XdmfUnstructuredGrid *grid, vtkUnstructuredGrid *dataSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Populates the VTK data set with the contents of the Xdmf grid.
static void CopyShape(XdmfUnstructuredGrid *grid, vtkUnstructuredGrid *dataSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Helper that does topology for XdmfToVTK.
static bool VTKToXdmfArray(vtkDataArray *vArray, XdmfArray *xArray, unsigned int rank=0, unsigned int *dims=nullptr)
Populates and Xdmf array corresponding to the VTK array it is given.
static void XdmfToVTKAttributes(XdmfSet *grid, vtkDataObject *dObject, vtkXdmf3ArrayKeeper *keeper=nullptr)
Populates the given VTK DataObject's attribute arrays with the selected arrays from the Xdmf Grid.
static void CopyShape(XdmfRectilinearGrid *grid, vtkRectilinearGrid *dataSet, vtkXdmf3ArrayKeeper *keeper=nullptr)
Helper that does topology for XdmfToVTK.
@ time
Definition: vtkX3D.h:497
@ name
Definition: vtkX3D.h:219
@ string
Definition: vtkX3D.h:490