VTK
DataSetConverters.h
Go to the documentation of this file.
1 //=============================================================================
2 //
3 // Copyright (c) Kitware, Inc.
4 // All rights reserved.
5 // See LICENSE.txt for details.
6 //
7 // This software is distributed WITHOUT ANY WARRANTY; without even
8 // the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
9 // PURPOSE. See the above copyright notice for more information.
10 //
11 // Copyright 2012 Sandia Corporation.
12 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
13 // the U.S. Government retains certain rights in this software.
14 //
15 //=============================================================================
16 
17 #ifndef daxToVtk_DataSetConverter_h
18 #define daxToVtk_DataSetConverter_h
19 
20 class vtkLine;
21 class vtkHexahedron;
22 class vtkQuad;
23 class vtkTetra;
24 class vtkTriangle;
25 class vtkVoxel;
26 class vtkWedge;
27 
28 #include "vtkCellData.h"
29 #include "vtkDataArray.h"
30 #include "vtkDataSet.h"
31 #include "vtkNew.h"
32 #include "vtkPointData.h"
33 #include "vtkPolyData.h"
34 #include "vtkUnstructuredGrid.h"
35 
36 #include <dax/cont/UnstructuredGrid.h>
37 #include <dax/cont/UniformGrid.h>
38 #include <dax/cont/ArrayHandle.h>
39 
40 #include "CellTypeToType.h"
41 #include <algorithm>
42 
43 
44 namespace daxToVtk
45 {
46 
47 namespace detail
48 {
49 //------------------------------------------------------------------------------
50 template<typename CellType>
52 {
53  //in no place do we in dax write the number of points are in the cell
54  //we don't want to that in the allocator. If the allocator does it
55  //we create an affinity between the thread the allocator is on and
56  //the memory, which will cause performance issues when we are in openMP
57 
58  //so instead we do it once we pull back to vtk
59  vtkIdType* raw_ids = cell->GetPointer();
60 
61  for(vtkIdType i=0; i < cell->GetNumberOfCells(); ++i)
62  {
63  raw_ids[i*(CellType::NUM_POINTS+1)]=CellType::NUM_POINTS;
64  }
65 }
66 
67 //------------------------------------------------------------------------------
68 template<typename CellType>
69 void setCells(vtkCellArray* cells, vtkPolyData* output)
70 {
71  //get the vtk cell type we are extracting
72  const VTKCellType cell_type = static_cast<VTKCellType>(CellType::VTKCellType);
73  if(cell_type == VTK_VERTEX)
74  {
75  output->SetVerts(cells);
76  }
77  else if(cell_type == VTK_LINE)
78  {
79  output->SetLines(cells);
80  }
81  else if(cell_type == VTK_TRIANGLE ||
82  cell_type == VTK_QUAD )
83  {
84  output->SetPolys(cells);
85  }
86 }
87 
88 //------------------------------------------------------------------------------
89 template<typename CellType>
91 {
92  //get the vtk cell type we are extracting
93  const VTKCellType cell_type = static_cast<VTKCellType>(CellType::VTKCellType);
94  output->SetCells(cell_type,cells);
95 }
96 
97 
98 //------------------------------------------------------------------------------
99 template<typename ContainerTag, typename GridType, typename OutputType>
100 void convertCells(ContainerTag, GridType& grid, OutputType* output)
101 {
102  //we are dealing with a container type whose memory wasn't allocated by
103  //vtk so we have to copy the data into a new vtk memory location just
104  //to be safe.
105  typedef typename ::daxToVtk::CellTypeToType<
106  typename GridType::CellTag > CellType;
107 
108  //determine amount of memory to allocate
109  const vtkIdType num_cells = grid.GetNumberOfCells();
110  const vtkIdType alloc_size = grid.GetCellConnections().GetNumberOfValues();
111 
112  //get the portal from the grid.
113  typedef typename GridType::CellConnectionsType::PortalConstControl
114  DaxPortalType;
115  DaxPortalType daxPortal = grid.GetCellConnections().GetPortalConstControl();
116 
117  //ask the vtkToDax allocator to make us memory
119  vtkCellArray* cells = alloc.allocate(num_cells+alloc_size);
120 
121  vtkIdType* cellPointer = cells->GetPointer();
122  vtkIdType index = 0;
123  for(vtkIdType i=0; i < num_cells; ++i)
124  {
125  *cellPointer = CellType::NUM_POINTS;
126  ++cellPointer;
127  for(vtkIdType j=0; j < CellType::NUM_POINTS; ++j, ++cellPointer, ++index)
128  {
129  *cellPointer = daxPortal.Get(index);
130  }
131  }
132 
133  daxToVtk::detail::setCells<CellType>(cells,output);
134 }
135 
136 //------------------------------------------------------------------------------
137 template<typename CellType, typename GridType, typename OutputType>
139  GridType& grid,
140  OutputType* output)
141 {
142  //in this use case the cell container is of vtk type so we
143  //can directly hook in and use the memory the container allocated
144  //for the output. This is really nice when working with TBB and OpenMP
145  //device adapters.
146  vtkCellArray* cells = grid.GetCellConnections().GetPortalControl().GetVtkData();
147 
148 
149  //to properly set the cells back into vtk we have to make
150  //sure that for each cell we will fill in the part which states
151  //how many points are in that cells
152  daxToVtk::detail::writeCellTags<CellType>(cells);
153 
154  daxToVtk::detail::setCells<CellType>(cells,output);
155 }
156 
157 //------------------------------------------------------------------------------
158 template<typename ContainerTag, typename GridType, typename OutputType>
159 void convertPoints(ContainerTag, GridType& grid, OutputType* output)
160 {
161  //we are dealing with a container type whose memory wasn't allocated by
162  //vtk so we have to copy the data into a new vtk memory location just
163  //to be safe.
164 
165  //determine amount of memory to allocate
166  const vtkIdType num_points = grid.GetNumberOfPoints();
167 
168  //ask vtkToDax to allocate the vtkPoints so it gets the float vs double
169  //settings correct
171  vtkPoints* points = alloc.allocate(num_points);
172 
173  dax::Vector3 *raw_pts = reinterpret_cast<dax::Vector3*>(
174  points->GetData()->GetVoidPointer(0));
175 
176  //get the coord portal from the grid.
177  typedef typename GridType::PointCoordinatesType::PortalConstControl
178  DaxPortalType;
179  DaxPortalType daxPortal = grid.GetPointCoordinates().GetPortalConstControl();
180 
181  std::copy(daxPortal.GetIteratorBegin(),
182  daxPortal.GetIteratorEnd(),
183  raw_pts);
184 
185  output->SetPoints( points );
186 }
187 
188 
189 //------------------------------------------------------------------------------
190 template<typename GridType, typename OutputType>
192  GridType& grid,
193  OutputType* output)
194 {
195  vtkPoints *p = grid.GetPointCoordinates().GetPortalControl().GetVtkData();
196  output->SetPoints(p);
197 }
198 
199 
200 } //namespace detail
201 
202 
203 //------------------------------------------------------------------------------
204 //convert a UniformGrid to vtkImageData
205 inline void dataSetConverter(const dax::cont::UniformGrid<>& grid,
206  vtkImageData* output)
207 {
208  dax::Vector3 origin = grid.GetOrigin();
209  dax::Vector3 spacing = grid.GetSpacing();
210  dax::Extent3 extent = grid.GetExtent();
211 
212  output->SetOrigin(origin[0],origin[1],origin[2]);
213  output->SetSpacing(spacing[0],spacing[1],spacing[2]);
214  output->SetExtent(extent.Min[0],extent.Max[0],
215  extent.Min[1],extent.Max[1],
216  extent.Min[2],extent.Max[2]);
217 }
218 
219 //convert a UnstructuredGrid to vtkUnstructuredGrid
220 template<typename CellType, typename TopoTag, typename PointTag>
221 inline void dataSetConverter(dax::cont::UnstructuredGrid<CellType,TopoTag,PointTag>& grid,
222  vtkUnstructuredGrid* output)
223 {
224  daxToVtk::detail::convertCells(TopoTag(),grid,output);
225  daxToVtk::detail::convertPoints(PointTag(),grid,output);
226 }
227 
228 //convert a UnstructuredGrid to vtkPolyData
229 template<typename CellType, typename TopoTag, typename PointTag>
230 inline void dataSetConverter(
231  dax::cont::UnstructuredGrid<CellType,TopoTag,PointTag>& grid,
232  vtkPolyData* output)
233 {
234  daxToVtk::detail::convertCells(TopoTag(),grid,output);
235  daxToVtk::detail::convertPoints(PointTag(),grid,output);
236 }
237 
238 template<typename FieldType>
239 void addCellData(vtkDataSet* output,
240  FieldType& outputArray,
241  const std::string& name)
242 {
243  vtkDataArray *data = outputArray.GetPortalControl().GetVtkData();
244  data->SetName(name.c_str());
245  output->GetCellData()->AddArray(data);
246 }
247 
248 
249 template<typename FieldType>
251  FieldType& outputArray,
252  const std::string& name)
253 {
254  vtkDataArray *data = outputArray.GetPortalControl().GetVtkData();
255  data->SetName(name.c_str());
256  vtkPointData *pd = output->GetPointData();
257  pd->AddArray(data);
258 }
259 
260 
261 }
262 #endif // daxToVtk_DataSetConverter_h
vtkIdType * GetPointer()
Get pointer to array of cell data.
Definition: vtkCellArray.h:227
pointer allocate(size_type n, self::const_pointer hint=0)
Definition: Allocators.h:53
int AddArray(vtkAbstractArray *array)
Add an array to the array list.
void dataSetConverter(const dax::cont::UniformGrid<> &grid, vtkImageData *output)
represent and manipulate point attribute data
Definition: vtkPointData.h:37
virtual void * GetVoidPointer(vtkIdType valueIdx)=0
Return a void pointer.
void SetLines(vtkCellArray *l)
Set the cell array defining lines.
void SetPolys(vtkCellArray *p)
Set the cell array defining polygons.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
vtkCellData * GetCellData()
Return a pointer to this dataset's cell data.
Definition: vtkDataSet.h:250
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:41
int vtkIdType
Definition: vtkType.h:287
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:85
virtual vtkIdType GetNumberOfCells()
Get the number of cells in the array.
VTKCellType
Definition: vtkCellType.h:43
void convertCells(ContainerTag, GridType &grid, OutputType *output)
virtual void SetName(const char *)
Set/get array's name.
vtkPointData * GetPointData()
Return a pointer to this dataset's point data.
Definition: vtkDataSet.h:256
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:47
cell represents a 1D line
Definition: vtkLine.h:35
virtual void SetExtent(int extent[6])
Set/Get the extent.
a cell that represents a 3D orthogonal parallelepiped
Definition: vtkVoxel.h:44
void writeCellTags(vtkCellArray *cell)
void convertPoints(ContainerTag, GridType &grid, OutputType *output)
topologically and geometrically regular array of data
Definition: vtkImageData.h:45
dataset represents arbitrary combinations of all possible cell types
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
void addPointData(vtkDataSet *output, FieldType &outputArray, const std::string &name)
void SetVerts(vtkCellArray *v)
Set the cell array defining vertices.
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:47
virtual void SetSpacing(double, double, double)
Set the spacing (width,height,length) of the cubical cells that compose the data set.
void SetCells(int type, vtkCellArray *cells)
Special methods specific to vtkUnstructuredGrid for defining the cells composing the dataset...
object to represent cell connectivity
Definition: vtkCellArray.h:50
a cell that represents a triangle
Definition: vtkTriangle.h:41
virtual void SetOrigin(double, double, double)
Set/Get the origin of the dataset.
void addCellData(vtkDataSet *output, FieldType &outputArray, const std::string &name)
vtkDataArray * GetData()
Definition: vtkPoints.h:69
void setCells(vtkCellArray *cells, vtkPolyData *output)
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:49
represent and manipulate 3D points
Definition: vtkPoints.h:39