Skip to content
Snippets Groups Projects
Commit c44f39f9 authored by Sujin Philip's avatar Sujin Philip
Browse files

Improved VTK-m to VTK Conversions

parent d5861562
No related branches found
No related tags found
No related merge requests found
......@@ -24,19 +24,25 @@
#include "vtkmDataArray.h"
// datasets we support
#include "vtkCellArray.h"
#include "vtkCellData.h"
#include "vtkCellTypes.h"
#include "vtkDataObject.h"
#include "vtkDataObjectTypes.h"
#include "vtkDataSetAttributes.h"
#include "vtkImageData.h"
#include "vtkNew.h"
#include "vtkPointData.h"
#include "vtkPolyData.h"
#include "vtkRectilinearGrid.h"
#include "vtkSmartPointer.h"
#include "vtkStructuredGrid.h"
#include "vtkUnstructuredGrid.h"
#include <vtkm/cont/ArrayHandle.h>
#include <vtkm/cont/ArrayHandleCartesianProduct.h>
#include <vtkm/cont/ArrayHandleUniformPointCoordinates.h>
#include <vtkm/cont/CellSetStructured.h>
#include <vtkm/cont/CoordinateSystem.hxx>
#include <vtkm/cont/Field.h>
......@@ -163,3 +169,142 @@ vtkm::cont::DataSet Convert(vtkDataSet* input, FieldsFlag fields)
}
} // namespace tovtkm
namespace fromvtkm
{
namespace
{
struct ComputeExtents
{
template <vtkm::IdComponent Dim>
void operator()(const vtkm::cont::CellSetStructured<Dim>& cs,
const vtkm::Id3& structuredCoordsDims, int extent[6]) const
{
auto extStart = cs.GetGlobalPointIndexStart();
for (int i = 0, ii = 0; i < 3; ++i)
{
if (structuredCoordsDims[i] > 1)
{
extent[2 * i] = vtkm::VecTraits<decltype(extStart)>::GetComponent(extStart, ii++);
extent[(2 * i) + 1] = extent[2 * i] + structuredCoordsDims[i] - 1;
}
else
{
extent[2 * i] = extent[(2 * i) + 1] = 0;
}
}
}
template <vtkm::IdComponent Dim>
void operator()(const vtkm::cont::CellSetStructured<Dim>& cs, int extent[6]) const
{
auto extStart = cs.GetGlobalPointIndexStart();
auto csDim = cs.GetPointDimensions();
for (int i = 0; i < Dim; ++i)
{
extent[2 * i] = vtkm::VecTraits<decltype(extStart)>::GetComponent(extStart, i);
extent[(2 * i) + 1] =
extent[2 * i] + vtkm::VecTraits<decltype(csDim)>::GetComponent(csDim, i) - 1;
}
for (int i = Dim; i < 3; ++i)
{
extent[2 * i] = extent[(2 * i) + 1] = 0;
}
}
};
} // anonymous namespace
void PassAttributesInformation(vtkDataSetAttributes* input, vtkDataSetAttributes* output)
{
for (int attribType = 0; attribType < vtkDataSetAttributes::NUM_ATTRIBUTES; attribType++)
{
vtkDataArray* attribute = input->GetAttribute(attribType);
if (attribute == nullptr)
{
continue;
}
output->SetActiveAttribute(attribute->GetName(), attribType);
}
}
bool Convert(const vtkm::cont::DataSet& vtkmOut, vtkRectilinearGrid* output, vtkDataSet* input)
{
using ListCellSetStructured = vtkm::List<vtkm::cont::CellSetStructured<1>,
vtkm::cont::CellSetStructured<2>, vtkm::cont::CellSetStructured<3> >;
auto cellSet = vtkmOut.GetCellSet().ResetCellSetList(ListCellSetStructured{});
using coordType =
vtkm::cont::ArrayHandleCartesianProduct<vtkm::cont::ArrayHandle<vtkm::FloatDefault>,
vtkm::cont::ArrayHandle<vtkm::FloatDefault>, vtkm::cont::ArrayHandle<vtkm::FloatDefault> >;
auto coordsArray = vtkm::cont::Cast<coordType>(vtkmOut.GetCoordinateSystem().GetData());
vtkSmartPointer<vtkDataArray> xArray =
Convert(vtkm::cont::make_FieldPoint("xArray", coordsArray.GetStorage().GetFirstArray()));
vtkSmartPointer<vtkDataArray> yArray =
Convert(vtkm::cont::make_FieldPoint("yArray", coordsArray.GetStorage().GetSecondArray()));
vtkSmartPointer<vtkDataArray> zArray =
Convert(vtkm::cont::make_FieldPoint("zArray", coordsArray.GetStorage().GetThirdArray()));
if (!xArray || !yArray || !zArray)
{
return false;
}
vtkm::Id3 dims(
xArray->GetNumberOfValues(), yArray->GetNumberOfValues(), zArray->GetNumberOfValues());
int extents[6];
vtkm::cont::CastAndCall(cellSet, ComputeExtents{}, dims, extents);
output->SetExtent(extents);
output->SetXCoordinates(xArray);
output->SetYCoordinates(yArray);
output->SetZCoordinates(zArray);
// Next we need to convert any extra fields from vtkm over to vtk
if (!fromvtkm::ConvertArrays(vtkmOut, output))
{
return false;
}
// Pass information about attributes.
PassAttributesInformation(input->GetPointData(), output->GetPointData());
PassAttributesInformation(input->GetCellData(), output->GetCellData());
return true;
}
bool Convert(const vtkm::cont::DataSet& vtkmOut, vtkStructuredGrid* output, vtkDataSet* input)
{
using ListCellSetStructured = vtkm::List<vtkm::cont::CellSetStructured<1>,
vtkm::cont::CellSetStructured<2>, vtkm::cont::CellSetStructured<3> >;
auto cellSet = vtkmOut.GetCellSet().ResetCellSetList(ListCellSetStructured{});
int extents[6];
vtkm::cont::CastAndCall(cellSet, ComputeExtents{}, extents);
vtkSmartPointer<vtkPoints> points = Convert(vtkmOut.GetCoordinateSystem());
if (!points)
{
return false;
}
output->SetExtent(extents);
output->SetPoints(points);
// Next we need to convert any extra fields from vtkm over to vtk
if (!fromvtkm::ConvertArrays(vtkmOut, output))
{
return false;
}
// Pass information about attributes.
PassAttributesInformation(input->GetPointData(), output->GetPointData());
PassAttributesInformation(input->GetCellData(), output->GetCellData());
return true;
}
} // namespace fromvtkm
......@@ -25,10 +25,12 @@
#include <vtkm/cont/DataSet.h>
class vtkDataSet;
class vtkDataSetAttributes;
class vtkImageData;
class vtkStructuredGrid;
class vtkPoints;
class vtkDataSet;
class vtkRectilinearGrid;
class vtkStructuredGrid;
namespace tovtkm
{
......@@ -46,4 +48,18 @@ VTKACCELERATORSVTKM_EXPORT
vtkm::cont::DataSet Convert(vtkDataSet* input, FieldsFlag fields = FieldsFlag::None);
}
namespace fromvtkm
{
VTKACCELERATORSVTKM_EXPORT
void PassAttributesInformation(vtkDataSetAttributes* input, vtkDataSetAttributes* output);
VTKACCELERATORSVTKM_EXPORT
bool Convert(const vtkm::cont::DataSet& vtkmOut, vtkRectilinearGrid* output, vtkDataSet* input);
VTKACCELERATORSVTKM_EXPORT
bool Convert(const vtkm::cont::DataSet& vtkmOut, vtkStructuredGrid* output, vtkDataSet* input);
}
#endif // vtkmlib_DataSetConverters_h
......@@ -18,12 +18,59 @@
#include <vtkm/cont/DataSetBuilderUniform.h>
#include "ArrayConverters.h"
#include "DataSetConverters.h"
#include "vtkCellData.h"
#include "vtkDataArray.h"
#include "vtkDataSetAttributes.h"
#include "vtkImageData.h"
#include "vtkPointData.h"
namespace
{
struct ComputeExtents
{
template <vtkm::IdComponent Dim>
void operator()(const vtkm::cont::CellSetStructured<Dim>& cs,
const vtkm::Id3& structuredCoordsDims, int extent[6]) const
{
auto extStart = cs.GetGlobalPointIndexStart();
for (int i = 0, ii = 0; i < 3; ++i)
{
if (structuredCoordsDims[i] > 1)
{
extent[2 * i] = vtkm::VecTraits<decltype(extStart)>::GetComponent(extStart, ii++);
extent[(2 * i) + 1] = extent[2 * i] + structuredCoordsDims[i] - 1;
}
else
{
extent[2 * i] = extent[(2 * i) + 1] = 0;
}
}
}
};
struct SetGlobalPointIndexStart
{
template <vtkm::IdComponent Dim, typename DynamicCellSetType>
void operator()(const vtkm::cont::CellSetStructured<Dim>&, const vtkm::Id3& structuredCoordsDims,
const int extent[6], DynamicCellSetType& dcs) const
{
typename vtkm::cont::CellSetStructured<Dim>::SchedulingRangeType extStart{};
for (int i = 0, ii = 0; i < 3; ++i)
{
if (structuredCoordsDims[i] > 1)
{
vtkm::VecTraits<decltype(extStart)>::SetComponent(extStart, ii++, extent[2 * i]);
}
}
vtkm::cont::Cast<vtkm::cont::CellSetStructured<Dim> >(dcs).SetGlobalPointIndexStart(extStart);
}
};
} // anonymous namespace
namespace tovtkm
{
......@@ -50,6 +97,11 @@ vtkm::cont::DataSet Convert(vtkImageData* input, FieldsFlag fields)
vtkm::cont::DataSet dataset = vtkm::cont::DataSetBuilderUniform::Create(dims, origin, spacing);
using ListCellSetStructured = vtkm::List<vtkm::cont::CellSetStructured<1>,
vtkm::cont::CellSetStructured<2>, vtkm::cont::CellSetStructured<3> >;
auto cellSet = dataset.GetCellSet().ResetCellSetList(ListCellSetStructured{});
vtkm::cont::CastAndCall(cellSet, SetGlobalPointIndexStart{}, dims, extent, dataset.GetCellSet());
ProcessFields(input, dataset, fields);
return dataset;
......@@ -90,15 +142,8 @@ bool Convert(
bool arraysConverted = fromvtkm::ConvertArrays(voutput, output);
// Pass information about attributes.
for (int attributeType = 0; attributeType < vtkDataSetAttributes::NUM_ATTRIBUTES; attributeType++)
{
vtkDataArray* attribute = input->GetPointData()->GetAttribute(attributeType);
if (attribute == nullptr)
{
continue;
}
output->GetPointData()->SetActiveAttribute(attribute->GetName(), attributeType);
}
PassAttributesInformation(input->GetPointData(), output->GetPointData());
PassAttributesInformation(input->GetCellData(), output->GetCellData());
return arraysConverted;
}
......@@ -115,8 +160,12 @@ bool Convert(const vtkm::cont::DataSet& voutput, vtkImageData* output, vtkDataSe
auto portal = points.GetPortalConstControl();
auto dim = portal.GetDimensions();
int extents[6] = { 0, static_cast<int>(dim[0] - 1), 0, static_cast<int>(dim[1] - 1), 0,
static_cast<int>(dim[2] - 1) };
int extents[6];
using ListCellSetStructured = vtkm::List<vtkm::cont::CellSetStructured<1>,
vtkm::cont::CellSetStructured<2>, vtkm::cont::CellSetStructured<3> >;
auto cellSet = voutput.GetCellSet().ResetCellSetList(ListCellSetStructured{});
vtkm::cont::CastAndCall(cellSet, ComputeExtents{}, dim, extents);
return Convert(voutput, extents, output, input);
}
......
......@@ -21,6 +21,7 @@
// datasets we support
#include "vtkCellArray.h"
#include "vtkCellData.h"
#include "vtkCellTypes.h"
#include "vtkDataObject.h"
#include "vtkDataObjectTypes.h"
......@@ -171,15 +172,8 @@ bool Convert(const vtkm::cont::DataSet& voutput, vtkPolyData* output, vtkDataSet
bool arraysConverted = ConvertArrays(voutput, output);
// Pass information about attributes.
for (int attributeType = 0; attributeType < vtkDataSetAttributes::NUM_ATTRIBUTES; attributeType++)
{
vtkDataArray* attribute = input->GetPointData()->GetAttribute(attributeType);
if (attribute == nullptr)
{
continue;
}
output->GetPointData()->SetActiveAttribute(attribute->GetName(), attributeType);
}
PassAttributesInformation(input->GetPointData(), output->GetPointData());
PassAttributesInformation(input->GetCellData(), output->GetCellData());
return arraysConverted;
}
......
......@@ -22,6 +22,7 @@
// datasets we support
#include "vtkCellArray.h"
#include "vtkCellData.h"
#include "vtkCellTypes.h"
#include "vtkDataObject.h"
#include "vtkDataObjectTypes.h"
......@@ -112,15 +113,8 @@ bool Convert(const vtkm::cont::DataSet& voutput, vtkUnstructuredGrid* output, vt
const bool arraysConverted = fromvtkm::ConvertArrays(voutput, output);
// Pass information about attributes.
for (int attributeType = 0; attributeType < vtkDataSetAttributes::NUM_ATTRIBUTES; attributeType++)
{
vtkDataArray* attribute = input->GetPointData()->GetAttribute(attributeType);
if (attribute == nullptr)
{
continue;
}
output->GetPointData()->SetActiveAttribute(attribute->GetName(), attributeType);
}
PassAttributesInformation(input->GetPointData(), output->GetPointData());
PassAttributesInformation(input->GetCellData(), output->GetCellData());
return arraysConverted;
}
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment