1 /*=========================================================================
2 
3   Program:   Visualization Toolkit
4   Module:    ImageDataConverter.cxx
5 
6   Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7   All rights reserved.
8   See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10      This software is distributed WITHOUT ANY WARRANTY; without even
11      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12      PURPOSE.  See the above copyright notice for more information.
13 
14 =========================================================================*/
15 #include "ImageDataConverter.h"
16 
17 #include <vtkm/cont/ArrayHandleUniformPointCoordinates.h>
18 #include <vtkm/cont/DataSetBuilderUniform.h>
19 
20 #include "ArrayConverters.h"
21 #include "DataSetConverters.h"
22 
23 #include "vtkCellData.h"
24 #include "vtkDataArray.h"
25 #include "vtkDataSetAttributes.h"
26 #include "vtkImageData.h"
27 #include "vtkPointData.h"
28 
29 namespace
30 {
31 
32 struct ComputeExtents
33 {
34   template <vtkm::IdComponent Dim>
operator ()__anon839d3e370111::ComputeExtents35   void operator()(const vtkm::cont::CellSetStructured<Dim>& cs,
36     const vtkm::Id3& structuredCoordsDims, int extent[6]) const
37   {
38     auto extStart = cs.GetGlobalPointIndexStart();
39     for (int i = 0, ii = 0; i < 3; ++i)
40     {
41       if (structuredCoordsDims[i] > 1)
42       {
43         extent[2 * i] = vtkm::VecTraits<decltype(extStart)>::GetComponent(extStart, ii++);
44         extent[(2 * i) + 1] = extent[2 * i] + structuredCoordsDims[i] - 1;
45       }
46       else
47       {
48         extent[2 * i] = extent[(2 * i) + 1] = 0;
49       }
50     }
51   }
52 };
53 
54 struct SetGlobalPointIndexStart
55 {
56   template <vtkm::IdComponent Dim, typename DynamicCellSetType>
operator ()__anon839d3e370111::SetGlobalPointIndexStart57   void operator()(const vtkm::cont::CellSetStructured<Dim>&, const vtkm::Id3& structuredCoordsDims,
58     const int extent[6], DynamicCellSetType& dcs) const
59   {
60     typename vtkm::cont::CellSetStructured<Dim>::SchedulingRangeType extStart{};
61     for (int i = 0, ii = 0; i < 3; ++i)
62     {
63       if (structuredCoordsDims[i] > 1)
64       {
65         vtkm::VecTraits<decltype(extStart)>::SetComponent(extStart, ii++, extent[2 * i]);
66       }
67     }
68     vtkm::cont::Cast<vtkm::cont::CellSetStructured<Dim>>(dcs).SetGlobalPointIndexStart(extStart);
69   }
70 };
71 
72 } // anonymous namespace
73 
74 namespace tovtkm
75 {
76 
77 //------------------------------------------------------------------------------
78 // convert an image data type
Convert(vtkImageData * input,FieldsFlag fields)79 vtkm::cont::DataSet Convert(vtkImageData* input, FieldsFlag fields)
80 {
81   int extent[6];
82   input->GetExtent(extent);
83   double vorigin[3];
84   input->GetOrigin(vorigin);
85   double vspacing[3];
86   input->GetSpacing(vspacing);
87   int vdims[3];
88   input->GetDimensions(vdims);
89 
90   vtkm::Vec<vtkm::FloatDefault, 3> origin(
91     static_cast<vtkm::FloatDefault>((static_cast<double>(extent[0]) * vspacing[0]) + vorigin[0]),
92     static_cast<vtkm::FloatDefault>((static_cast<double>(extent[2]) * vspacing[1]) + vorigin[1]),
93     static_cast<vtkm::FloatDefault>((static_cast<double>(extent[4]) * vspacing[2]) + vorigin[2]));
94   vtkm::Vec<vtkm::FloatDefault, 3> spacing(static_cast<vtkm::FloatDefault>(vspacing[0]),
95     static_cast<vtkm::FloatDefault>(vspacing[1]), static_cast<vtkm::FloatDefault>(vspacing[2]));
96   vtkm::Id3 dims(vdims[0], vdims[1], vdims[2]);
97 
98   vtkm::cont::DataSet dataset = vtkm::cont::DataSetBuilderUniform::Create(dims, origin, spacing);
99 
100   using ListCellSetStructured = vtkm::List<vtkm::cont::CellSetStructured<1>,
101     vtkm::cont::CellSetStructured<2>, vtkm::cont::CellSetStructured<3>>;
102   auto cellSet = dataset.GetCellSet().ResetCellSetList(ListCellSetStructured{});
103   vtkm::cont::CastAndCall(cellSet, SetGlobalPointIndexStart{}, dims, extent, dataset.GetCellSet());
104 
105   ProcessFields(input, dataset, fields);
106 
107   return dataset;
108 }
109 
110 } // tovtkm
111 
112 namespace fromvtkm
113 {
114 
Convert(const vtkm::cont::DataSet & voutput,int extents[6],vtkImageData * output,vtkDataSet * input)115 bool Convert(
116   const vtkm::cont::DataSet& voutput, int extents[6], vtkImageData* output, vtkDataSet* input)
117 {
118   vtkm::cont::CoordinateSystem const& cs = voutput.GetCoordinateSystem();
119   if (!cs.GetData().IsType<vtkm::cont::ArrayHandleUniformPointCoordinates>())
120   {
121     return false;
122   }
123 
124   auto points = cs.GetData().AsArrayHandle<vtkm::cont::ArrayHandleUniformPointCoordinates>();
125   auto portal = points.ReadPortal();
126 
127   auto origin = portal.GetOrigin();
128   auto spacing = portal.GetSpacing();
129   auto dim = portal.GetDimensions();
130   VTKM_ASSERT((extents[1] - extents[0] + 1) == dim[0] && (extents[3] - extents[2] + 1) == dim[1] &&
131     (extents[5] - extents[4] + 1) == dim[2]);
132 
133   origin[0] -= static_cast<vtkm::FloatDefault>(extents[0]) * spacing[0];
134   origin[1] -= static_cast<vtkm::FloatDefault>(extents[2]) * spacing[1];
135   origin[2] -= static_cast<vtkm::FloatDefault>(extents[4]) * spacing[2];
136 
137   output->SetExtent(extents);
138   output->SetOrigin(origin[0], origin[1], origin[2]);
139   output->SetSpacing(spacing[0], spacing[1], spacing[2]);
140 
141   // Next we need to convert any extra fields from vtkm over to vtk
142   bool arraysConverted = fromvtkm::ConvertArrays(voutput, output);
143 
144   // Pass information about attributes.
145   PassAttributesInformation(input->GetPointData(), output->GetPointData());
146   PassAttributesInformation(input->GetCellData(), output->GetCellData());
147 
148   return arraysConverted;
149 }
150 
Convert(const vtkm::cont::DataSet & voutput,vtkImageData * output,vtkDataSet * input)151 bool Convert(const vtkm::cont::DataSet& voutput, vtkImageData* output, vtkDataSet* input)
152 {
153   vtkm::cont::CoordinateSystem const& cs = voutput.GetCoordinateSystem();
154   if (!cs.GetData().IsType<vtkm::cont::ArrayHandleUniformPointCoordinates>())
155   {
156     return false;
157   }
158 
159   auto points = cs.GetData().AsArrayHandle<vtkm::cont::ArrayHandleUniformPointCoordinates>();
160   auto portal = points.ReadPortal();
161 
162   auto dim = portal.GetDimensions();
163   int extents[6];
164   using ListCellSetStructured = vtkm::List<vtkm::cont::CellSetStructured<1>,
165     vtkm::cont::CellSetStructured<2>, vtkm::cont::CellSetStructured<3>>;
166   auto cellSet = voutput.GetCellSet().ResetCellSetList(ListCellSetStructured{});
167   vtkm::cont::CastAndCall(cellSet, ComputeExtents{}, dim, extents);
168 
169   return Convert(voutput, extents, output, input);
170 }
171 
172 } // fromvtkm
173