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 vtkToDax_DataSetConverter_h
18#define vtkToDax_DataSetConverter_h
19
20//datasets we support
21#include "vtkDataObjectTypes.h"
22#include "vtkCellTypes.h"
23#include "vtkCellArray.h"
24#include "vtkImageData.h"
25#include "vtkUniformGrid.h"
26#include "vtkUnstructuredGrid.h"
27#include "vtkStructuredGrid.h"
28
29#include <dax/cont/ArrayHandle.h>
30#include <dax/cont/UniformGrid.h>
31#include <dax/cont/UnstructuredGrid.h>
32
33#include "CellTypeToType.h"
34#include "DataSetTypeToType.h"
35
36#include "Portals.h"
37#include "Containers.h"
38
39
40namespace vtkToDax
41{
42
43template<typename CellType>
45 std::vector<dax::Id>& topo)
46{
47 enum{NUM_POINTS_IN_CELL=CellType::NUM_POINTS};
49 topo.reserve(size*NUM_POINTS_IN_CELL); //reserver room so we don't have to realloc
50 vtkCellArray *cells = input->GetCells();
51
52 vtkIdType npts, *pts;
53 cells->InitTraversal();
54 while(cells->GetNextCell(npts,pts))
55 {
56 std::copy(pts,pts+npts,std::back_inserter(topo));
57 }
58}
59
60
61
62//convert an image data type
63template<typename VTKDataSetType>
64inline typename VTKDataSetType::DaxDataSetType dataSetConverter(
65 vtkImageData* input,
66 VTKDataSetType)
67 {
68 typedef typename VTKDataSetType::DaxDataSetType DataSet;
69 double origin[3];input->GetOrigin(origin);
70 double spacing[3];input->GetSpacing(spacing);
71 int extent[6];input->GetExtent(extent);
72
73 //this would be image data
74 DataSet output;
75 output.SetOrigin(dax::make_Vector3(origin[0],origin[1],origin[2]));
76 output.SetSpacing(dax::make_Vector3(spacing[0],spacing[1],spacing[2]));
77 output.SetExtent(dax::make_Id3(extent[0],extent[2],extent[4]),
78 dax::make_Id3(extent[1],extent[3],extent[5]));
79 return output;
80 }
81
82//convert an uniform grid type
83template<typename VTKDataSetType>
84inline typename VTKDataSetType::DaxDataSetType dataSetConverter(
85 vtkUniformGrid* input,
86 VTKDataSetType)
87 {
88 typedef typename VTKDataSetType::DaxDataSetType DataSet;
89 double origin[3];input->GetOrigin(origin);
90 double spacing[3];input->GetSpacing(spacing);
91 int extent[6];input->GetExtent(extent);
92
93 //this would be image data
94 DataSet output;
95 output.SetOrigin(dax::make_Vector3(origin[0],origin[1],origin[2]));
96 output.SetSpacing(dax::make_Vector3(spacing[0],spacing[1],spacing[2]));
97 output.SetExtent(dax::make_Id3(extent[0],extent[2],extent[4]),
98 dax::make_Id3(extent[1],extent[3],extent[5]));
99 return output;
100 }
101
102//convert an unstructured grid type
103template<typename VTKDataSetType>
104inline typename VTKDataSetType::DaxDataSetType dataSetConverter(
105 vtkUnstructuredGrid* input,
106 VTKDataSetType)
107 {
108 //we convert to a hexahedron unstruction grid
109 //this uses a vtkCellArrayContainer to get the needed topology information
110 typedef typename VTKDataSetType::DaxDataSetType DataSet;
111 typedef typename VTKDataSetType::CellTypeToType CellTypeToType;
112
113 static const int NUM_POINTS = VTKDataSetType::CellTypeToType::NUM_POINTS;
114
115
116 dax::cont::ArrayHandle<dax::Vector3,vtkToDax::vtkPointsContainerTag>
118 input->GetNumberOfPoints()));
119
120 //
121 dax::cont::ArrayHandle<dax::Id,vtkToDax::vtkTopologyContainerTag<CellTypeToType> >
123 input->GetNumberOfCells()*NUM_POINTS));
124
125 return DataSet(topoHandle,pointsHandle);
126 }
127}
128
129
130#endif // vtkToDax_DataSetConverter_h
object to represent cell connectivity
Definition: vtkCellArray.h:51
int GetNextCell(vtkIdType &npts, vtkIdType *&pts)
A cell traversal methods that is more efficient than vtkDataSet traversal methods.
Definition: vtkCellArray.h:357
void InitTraversal()
A cell traversal methods that is more efficient than vtkDataSet traversal methods.
Definition: vtkCellArray.h:102
topologically and geometrically regular array of data
Definition: vtkImageData.h:46
virtual double * GetOrigin()
virtual int * GetExtent()
virtual double * GetSpacing()
virtual vtkPoints * GetPoints()
vtkIdType GetNumberOfPoints() override
See vtkDataSet for additional information.
Definition: vtkPointSet.h:164
image data with blanking
dataset represents arbitrary combinations of all possible cell types
vtkCellArray * GetCells()
vtkIdType GetNumberOfCells() override
Determine the number of cells composing the dataset.
void convertTopology(vtkUnstructuredGrid *input, std::vector< dax::Id > &topo)
VTKDataSetType::DaxDataSetType dataSetConverter(vtkImageData *input, VTKDataSetType)
@ extent
Definition: vtkX3D.h:345
@ spacing
Definition: vtkX3D.h:481
@ size
Definition: vtkX3D.h:253
int vtkIdType
Definition: vtkType.h:287