VTK  9.1.0
vtkXdmf3DataSet.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkXdmf3DataSet.h
5 Language: C++
6
7 Copyright (c) 1993-2002 Ken Martin, Will Schroeder, Bill Lorensen
8 All rights reserved.
9 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
10
11 This software is distributed WITHOUT ANY WARRANTY; without even
12 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
13 PURPOSE. See the above copyright notice for more information.
14
15=========================================================================*/
27#ifndef vtkXdmf3DataSet_h
28#define vtkXdmf3DataSet_h
29
30#include "vtkIOXdmf3Module.h" // For export macro
31
32// clang-format off
33#include "vtk_xdmf3.h"
34#include VTKXDMF3_HEADER(core/XdmfSharedPtr.hpp)
35// clang-format on
36
37#include <string> //Needed only for XdmfArray::getName :(
38
41class XdmfArray;
42class XdmfAttribute;
43class vtkDataArray;
44class XdmfGrid;
45class vtkDataObject;
46class XdmfSet;
47class vtkDataSet;
48class XdmfTopologyType;
49class XdmfRegularGrid;
50class vtkImageData;
51class XdmfRectilinearGrid;
53class XdmfCurvilinearGrid;
55class XdmfUnstructuredGrid;
57class vtkPointSet;
58class XdmfGraph;
61class XdmfDomain;
62
63class VTKIOXDMF3_EXPORT vtkXdmf3DataSet
64{
65public:
66 // Common
67
71 static vtkDataArray* XdmfToVTKArray(XdmfArray* xArray,
72 std::string attrName, // TODO: needed because XdmfArray::getName() misbehaves
73 unsigned int preferredComponents = 0, vtkXdmf3ArrayKeeper* keeper = nullptr);
74
78 static bool VTKToXdmfArray(
79 vtkDataArray* vArray, XdmfArray* xArray, unsigned int rank = 0, unsigned int* dims = nullptr);
80
86 vtkXdmf3ArraySelection* cselection, vtkXdmf3ArraySelection* pselection, XdmfGrid* grid,
87 vtkDataObject* dObject, vtkXdmf3ArrayKeeper* keeper = nullptr);
88
93 static void VTKToXdmfAttributes(vtkDataObject* dObject, XdmfGrid* grid);
94
96
99 static unsigned int GetNumberOfPointsPerCell(int vtk_cell_type, bool& fail);
100 static int GetVTKCellType(shared_ptr<const XdmfTopologyType> topologyType);
101 static int GetXdmfCellType(int vtkType);
103
105
108 static void SetTime(XdmfGrid* grid, double hasTime, double time);
109 static void SetTime(XdmfGraph* graph, double hasTime, double time);
111
112 // vtkXdmf3RegularGrid
113
117 static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
118 vtkXdmf3ArraySelection* pselection, XdmfRegularGrid* grid, vtkImageData* dataSet,
119 vtkXdmf3ArrayKeeper* keeper = nullptr);
120
124 static void CopyShape(
125 XdmfRegularGrid* grid, vtkImageData* dataSet, vtkXdmf3ArrayKeeper* keeper = nullptr);
126
130 static void VTKToXdmf(
131 vtkImageData* dataSet, XdmfDomain* domain, bool hasTime, double time, const char* name = 0);
132
133 // vtkXdmf3RectilinearGrid
137 static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
138 vtkXdmf3ArraySelection* pselection, XdmfRectilinearGrid* grid, vtkRectilinearGrid* dataSet,
139 vtkXdmf3ArrayKeeper* keeper = nullptr);
140
144 static void CopyShape(
145 XdmfRectilinearGrid* grid, vtkRectilinearGrid* dataSet, vtkXdmf3ArrayKeeper* keeper = nullptr);
146
150 static void VTKToXdmf(vtkRectilinearGrid* dataSet, XdmfDomain* domain, bool hasTime, double time,
151 const char* name = 0);
152
153 // vtkXdmf3CurvilinearGrid
157 static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
158 vtkXdmf3ArraySelection* pselection, XdmfCurvilinearGrid* grid, vtkStructuredGrid* dataSet,
159 vtkXdmf3ArrayKeeper* keeper = nullptr);
160
164 static void CopyShape(
165 XdmfCurvilinearGrid* grid, vtkStructuredGrid* dataSet, vtkXdmf3ArrayKeeper* keeper = nullptr);
166
170 static void VTKToXdmf(vtkStructuredGrid* dataSet, XdmfDomain* domain, bool hasTime, double time,
171 const char* name = 0);
172
173 // vtkXdmf3UnstructuredGrid
177 static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
178 vtkXdmf3ArraySelection* pselection, XdmfUnstructuredGrid* grid, vtkUnstructuredGrid* dataSet,
179 vtkXdmf3ArrayKeeper* keeper = nullptr);
180
184 static void CopyShape(XdmfUnstructuredGrid* grid, vtkUnstructuredGrid* dataSet,
185 vtkXdmf3ArrayKeeper* keeper = nullptr);
186
190 static void VTKToXdmf(
191 vtkPointSet* dataSet, XdmfDomain* domain, bool hasTime, double time, const char* name = 0);
192
193 // vtkXdmf3Graph
197 static void XdmfToVTK(vtkXdmf3ArraySelection* fselection, vtkXdmf3ArraySelection* cselection,
198 vtkXdmf3ArraySelection* pselection, XdmfGraph* grid, vtkMutableDirectedGraph* dataSet,
199 vtkXdmf3ArrayKeeper* keeper = nullptr);
200
204 static void VTKToXdmf(
205 vtkDirectedGraph* dataSet, XdmfDomain* domain, bool hasTime, double time, const char* name = 0);
206
207 // Side Sets
208
214 /*
215 vtkXdmf3ArraySelection *fselection,
216 vtkXdmf3ArraySelection *cselection,
217 vtkXdmf3ArraySelection *pselection,
218 */
219 XdmfSet* grid, vtkDataObject* dObject, vtkXdmf3ArrayKeeper* keeper = nullptr);
220
225 static void XdmfSubsetToVTK(XdmfGrid* grid, unsigned int setnum, vtkDataSet* dataSet,
226 vtkUnstructuredGrid* subSet, vtkXdmf3ArrayKeeper* keeper = nullptr);
227
233 static int GetVTKFiniteElementCellType(unsigned int element_degree,
234 const std::string& element_family, shared_ptr<const XdmfTopologyType> topologyType);
235
248 shared_ptr<XdmfAttribute> xmfAttribute, vtkDataArray* array, XdmfGrid* grid,
249 vtkXdmf3ArrayKeeper* keeper = nullptr);
250};
251
252#endif
253// VTK-HeaderTest-Exclude: vtkXdmf3DataSet.h
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
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:157
An editable directed graph.
concrete class for storing a set of points
Definition: vtkPointSet.h:106
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 void VTKToXdmf(vtkRectilinearGrid *dataSet, XdmfDomain *domain, bool hasTime, double time, const char *name=0)
Populates the Xdmf Grid with the contents of the VTK data set.
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 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 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(vtkImageData *dataSet, XdmfDomain *domain, bool hasTime, double time, const char *name=0)
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=0)
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 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 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 VTKToXdmf(vtkDirectedGraph *dataSet, XdmfDomain *domain, bool hasTime, double time, const char *name=0)
Populates the Xdmf Grid with the contents of the VTK data set.
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 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 VTKToXdmf(vtkStructuredGrid *dataSet, XdmfDomain *domain, bool hasTime, double time, const char *name=0)
Populates the Xdmf Grid with the contents of the VTK data set.
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:503
@ name
Definition: vtkX3D.h:225
@ string
Definition: vtkX3D.h:496