Implement CellToPoint for CellSetPermutation

This commit is contained in:
Sujin Philip 2017-07-21 12:19:19 -04:00
parent 46aa74519a
commit 6ce8218d09
6 changed files with 513 additions and 101 deletions

@ -22,30 +22,157 @@
#include <vtkm/CellShape.h>
#include <vtkm/CellTraits.h>
#include <vtkm/cont/ArrayHandleCast.h>
#include <vtkm/cont/ArrayHandlePermutation.h>
#include <vtkm/cont/CellSet.h>
#include <vtkm/cont/internal/ConnectivityExplicitInternals.h>
#include <vtkm/internal/ConnectivityStructuredInternals.h>
#include <vtkm/worklet/DispatcherMapTopology.h>
#include <vtkm/worklet/WorkletMapTopology.h>
#include <vtkm/exec/ConnectivityPermuted.h>
#ifndef VTKM_DEFAULT_CELLSET_PERMUTATION_STORAGE_TAG
#define VTKM_DEFAULT_CELLSET_PERMUTATION_STORAGE_TAG VTKM_DEFAULT_STORAGE_TAG
#endif
namespace vtkm
{
namespace cont
{
template <typename OriginalCellSet, typename ValidCellArrayHandleType>
class CellSetPermutation;
namespace internal
{
template <typename OriginalCellSet, typename PermutationArrayHandleType>
class CellSetGeneralPermutation : public CellSet
struct WriteConnectivity : public vtkm::worklet::WorkletMapPointToCell
{
typedef void ControlSignature(CellSetIn cellset,
FieldInCell<IdType> offset,
WholeArrayOut<> connectivity);
typedef void ExecutionSignature(PointCount, PointIndices, _2, _3);
typedef _1 InputDomain;
template <typename PointIndicesType, typename OutPortalType>
VTKM_EXEC void operator()(vtkm::IdComponent pointcount,
const PointIndicesType& pointIndices,
vtkm::Id offset,
OutPortalType& connectivity) const
{
for (vtkm::IdComponent i = 0; i < pointcount; ++i)
{
connectivity.Set(offset++, pointIndices[i]);
}
}
};
// default for CellSetExplicit and CellSetSingleType
template <typename CellSetPermutationType>
class PointToCell
{
private:
using NumIndicesArrayType = decltype(make_ArrayHandlePermutation(
std::declval<CellSetPermutationType>().GetValidCellIds(),
std::declval<CellSetPermutationType>().GetFullCellSet().GetNumIndicesArray(
vtkm::TopologyElementTagPoint(),
vtkm::TopologyElementTagCell())));
public:
using ConnectivityArrays = vtkm::cont::internal::ConnectivityExplicitInternals<
VTKM_DEFAULT_STORAGE_TAG, // dummy, shapes array is not used
typename NumIndicesArrayType::StorageTag>;
template <typename Device>
static ConnectivityArrays Get(const CellSetPermutationType& cellset, Device)
{
ConnectivityArrays conn;
conn.NumIndices =
NumIndicesArrayType(cellset.GetValidCellIds(),
cellset.GetFullCellSet().GetNumIndicesArray(
vtkm::TopologyElementTagPoint(), vtkm::TopologyElementTagCell()));
vtkm::Id connectivityLength = vtkm::cont::DeviceAdapterAlgorithm<Device>::ScanExclusive(
vtkm::cont::make_ArrayHandleCast(conn.NumIndices, vtkm::Id()), conn.IndexOffsets);
conn.Connectivity.Allocate(connectivityLength);
vtkm::worklet::DispatcherMapTopology<WriteConnectivity, Device>().Invoke(
cellset, conn.IndexOffsets, conn.Connectivity);
return conn;
}
};
// specialization for CellSetStructured
template <vtkm::IdComponent DIMENSION, typename PermutationArrayHandleType>
class PointToCell<CellSetPermutation<CellSetStructured<DIMENSION>, PermutationArrayHandleType>>
{
private:
using CellSetPermutationType =
CellSetPermutation<CellSetStructured<DIMENSION>, PermutationArrayHandleType>;
public:
using ConnectivityArrays = vtkm::cont::internal::ConnectivityExplicitInternals<
VTKM_DEFAULT_STORAGE_TAG, // dummy, shapes array is not used
typename vtkm::cont::ArrayHandleConstant<vtkm::IdComponent>::StorageTag,
VTKM_DEFAULT_STORAGE_TAG,
typename vtkm::cont::ArrayHandleCounting<vtkm::Id>::StorageTag>;
template <typename Device>
static ConnectivityArrays Get(const CellSetPermutationType& cellset, Device)
{
vtkm::Id numberOfCells = cellset.GetNumberOfCells();
vtkm::IdComponent numPointsInCell =
vtkm::internal::ConnectivityStructuredInternals<DIMENSION>::NUM_POINTS_IN_CELL;
vtkm::Id connectivityLength = numberOfCells * numPointsInCell;
ConnectivityArrays conn;
conn.NumIndices = make_ArrayHandleConstant(numPointsInCell, numberOfCells);
conn.IndexOffsets = ArrayHandleCounting<vtkm::Id>(0, numPointsInCell, numberOfCells);
conn.Connectivity.Allocate(connectivityLength);
vtkm::worklet::DispatcherMapTopology<WriteConnectivity, Device>().Invoke(
cellset, conn.IndexOffsets, conn.Connectivity);
return conn;
}
};
template <typename CellSetPermutationType>
struct CellSetPermutationTraits;
template <typename OriginalCellSet_, typename PermutationArrayHandleType_>
struct CellSetPermutationTraits<CellSetPermutation<OriginalCellSet_, PermutationArrayHandleType_>>
{
using OriginalCellSet = OriginalCellSet_;
using PermutationArrayHandleType = PermutationArrayHandleType_;
};
template <typename OriginalCellSet_,
typename OriginalPermutationArrayHandleType,
typename PermutationArrayHandleType_>
struct CellSetPermutationTraits<
CellSetPermutation<CellSetPermutation<OriginalCellSet_, OriginalPermutationArrayHandleType>,
PermutationArrayHandleType_>>
{
using PreviousCellSet = CellSetPermutation<OriginalCellSet_, OriginalPermutationArrayHandleType>;
using PermutationArrayHandleType = vtkm::cont::ArrayHandlePermutation<
PermutationArrayHandleType_,
typename CellSetPermutationTraits<PreviousCellSet>::PermutationArrayHandleType>;
using OriginalCellSet = typename CellSetPermutationTraits<PreviousCellSet>::OriginalCellSet;
using Superclass = CellSetPermutation<OriginalCellSet, PermutationArrayHandleType>;
};
} // internal
template <typename OriginalCellSetType_,
typename PermutationArrayHandleType_ =
vtkm::cont::ArrayHandle<vtkm::Id, VTKM_DEFAULT_CELLSET_PERMUTATION_STORAGE_TAG>>
class CellSetPermutation : public CellSet
{
public:
using OriginalCellSetType = OriginalCellSetType_;
using PermutationArrayHandleType = PermutationArrayHandleType_;
VTKM_CONT
CellSetGeneralPermutation(const PermutationArrayHandleType& validCellIds,
const OriginalCellSet& cellset,
const std::string& name)
CellSetPermutation(const PermutationArrayHandleType& validCellIds,
const OriginalCellSetType& cellset,
const std::string& name = std::string())
: CellSet(name)
, ValidCellIds(validCellIds)
, FullCellSet(cellset)
@ -53,7 +180,7 @@ public:
}
VTKM_CONT
CellSetGeneralPermutation(const std::string& name)
CellSetPermutation(const std::string& name = std::string())
: CellSet(name)
, ValidCellIds()
, FullCellSet()
@ -61,130 +188,174 @@ public:
}
VTKM_CONT
vtkm::Id GetNumberOfCells() const VTKM_OVERRIDE { return this->ValidCellIds.GetNumberOfValues(); }
const OriginalCellSetType& GetFullCellSet() const { return this->FullCellSet; }
VTKM_CONT
vtkm::Id GetNumberOfPoints() const VTKM_OVERRIDE { return this->FullCellSet.GetNumberOfPoints(); }
const PermutationArrayHandleType& GetValidCellIds() const { return this->ValidCellIds; }
vtkm::Id GetNumberOfFaces() const VTKM_OVERRIDE { return -1; }
VTKM_CONT
vtkm::Id GetNumberOfCells() const override { return this->ValidCellIds.GetNumberOfValues(); }
vtkm::Id GetNumberOfEdges() const VTKM_OVERRIDE { return -1; }
VTKM_CONT
vtkm::Id GetNumberOfPoints() const override { return this->FullCellSet.GetNumberOfPoints(); }
VTKM_CONT
vtkm::Id GetNumberOfFaces() const override { return -1; }
VTKM_CONT
vtkm::Id GetNumberOfEdges() const override { return -1; }
VTKM_CONT
vtkm::IdComponent GetNumberOfPointsInCell(vtkm::Id cellIndex) const
{
return this->FullCellSet.GetNumberOfPointsInCell(
this->ValidCellIds.GetPortalConstControl().Get(cellIndex));
}
//This is the way you can fill the memory from another system without copying
VTKM_CONT
void Fill(const PermutationArrayHandleType& validCellIds, const OriginalCellSet& cellset)
void Fill(const PermutationArrayHandleType& validCellIds, const OriginalCellSetType& cellset)
{
ValidCellIds = validCellIds;
FullCellSet = cellset;
this->ValidCellIds = validCellIds;
this->FullCellSet = cellset;
}
template <typename TopologyElement>
VTKM_CONT vtkm::Id GetSchedulingRange(TopologyElement) const
VTKM_CONT vtkm::Id GetSchedulingRange(vtkm::TopologyElementTagCell) const
{
VTKM_IS_TOPOLOGY_ELEMENT_TAG(TopologyElement);
return this->ValidCellIds.GetNumberOfValues();
}
VTKM_CONT vtkm::Id GetSchedulingRange(vtkm::TopologyElementTagPoint) const
{
return this->FullCellSet.GetNumberOfPoints();
}
template <typename Device, typename FromTopology, typename ToTopology>
struct ExecutionTypes
struct ExecutionTypes;
template <typename Device>
struct ExecutionTypes<Device, vtkm::TopologyElementTagPoint, vtkm::TopologyElementTagCell>
{
VTKM_IS_DEVICE_ADAPTER_TAG(Device);
VTKM_IS_TOPOLOGY_ELEMENT_TAG(FromTopology);
VTKM_IS_TOPOLOGY_ELEMENT_TAG(ToTopology);
typedef typename PermutationArrayHandleType::template ExecutionTypes<Device>::PortalConst
ExecPortalType;
using ExecPortalType =
typename PermutationArrayHandleType::template ExecutionTypes<Device>::PortalConst;
using OrigExecObjectType = typename OriginalCellSetType::template ExecutionTypes<
Device,
vtkm::TopologyElementTagPoint,
vtkm::TopologyElementTagCell>::ExecObjectType;
typedef typename OriginalCellSet::template ExecutionTypes<Device, FromTopology, ToTopology>::
ExecObjectType OrigExecObjectType;
typedef vtkm::exec::ConnectivityPermuted<ExecPortalType, OrigExecObjectType> ExecObjectType;
using ExecObjectType =
vtkm::exec::ConnectivityPermutedPointToCell<ExecPortalType, OrigExecObjectType>;
};
template <typename Device, typename FromTopology, typename ToTopology>
typename ExecutionTypes<Device, FromTopology, ToTopology>::ExecObjectType
PrepareForInput(Device, FromTopology, ToTopology) const
template <typename Device>
struct ExecutionTypes<Device, vtkm::TopologyElementTagCell, vtkm::TopologyElementTagPoint>
{
// Developer's note: I looked into supporting cell to point connectivity in a permutation cell
// set and found it to be complex. It is not straightforward to implement this on top of the
// original cell set's cell to point because points could be attached to cells that do not
// exist in the permuted topology. Ultimately, you will probably have to rebuild these
// connections in a way very similar to how CellSetExplicit already does it. In fact, the
// easiest implementation will probably be to just convert to a CellSetExplicit and use that.
// In fact, it may be possible to change this whole implementation to just be a subclass of
// CellSetExplicit with some fancy arrays for its point to cell arrays.
throw vtkm::cont::ErrorBadType(
"CellSetPermutation currently only supports point to cell connectivity. "
"To support other connectivity, convert to an explicit grid with the CellDeepCopy "
"worklet or the CleanGrid filter.");
VTKM_IS_DEVICE_ADAPTER_TAG(Device);
using ConnectiviyPortalType =
typename vtkm::cont::ArrayHandle<vtkm::Id>::template ExecutionTypes<Device>::PortalConst;
using NumIndicesPortalType = typename vtkm::cont::ArrayHandle<
vtkm::IdComponent>::template ExecutionTypes<Device>::PortalConst;
using IndexOffsetPortalType =
typename vtkm::cont::ArrayHandle<vtkm::Id>::template ExecutionTypes<Device>::PortalConst;
using ExecObjectType = vtkm::exec::ConnectivityPermutedCellToPoint<ConnectiviyPortalType,
NumIndicesPortalType,
IndexOffsetPortalType>;
};
template <typename Device>
VTKM_CONT typename ExecutionTypes<Device,
vtkm::TopologyElementTagPoint,
vtkm::TopologyElementTagCell>::ExecObjectType
PrepareForInput(Device device,
vtkm::TopologyElementTagPoint from,
vtkm::TopologyElementTagCell to) const
{
using ConnectivityType = typename ExecutionTypes<Device,
vtkm::TopologyElementTagPoint,
vtkm::TopologyElementTagCell>::ExecObjectType;
return ConnectivityType(this->ValidCellIds.PrepareForInput(device),
this->FullCellSet.PrepareForInput(device, from, to));
}
template <typename Device>
typename ExecutionTypes<Device,
vtkm::TopologyElementTagPoint,
vtkm::TopologyElementTagCell>::ExecObjectType
PrepareForInput(Device d, vtkm::TopologyElementTagPoint f, vtkm::TopologyElementTagCell t) const
VTKM_CONT typename ExecutionTypes<Device,
vtkm::TopologyElementTagCell,
vtkm::TopologyElementTagPoint>::ExecObjectType
PrepareForInput(Device device, vtkm::TopologyElementTagCell, vtkm::TopologyElementTagPoint) const
{
using FromTopology = vtkm::TopologyElementTagPoint;
using ToTopology = vtkm::TopologyElementTagCell;
using ConnectivityType =
typename ExecutionTypes<Device, FromTopology, ToTopology>::ExecObjectType;
return ConnectivityType(this->ValidCellIds.PrepareForInput(d),
this->FullCellSet.PrepareForInput(d, f, t));
if (!this->CellToPoint.ElementsValid)
{
auto pointToCell = internal::PointToCell<CellSetPermutation>::Get(*this, device);
internal::ComputeCellToPointConnectivity(
this->CellToPoint, pointToCell, this->GetNumberOfPoints(), device);
this->CellToPoint.BuildIndexOffsets(device);
}
using ConnectivityType = typename ExecutionTypes<Device,
vtkm::TopologyElementTagCell,
vtkm::TopologyElementTagPoint>::ExecObjectType;
return ConnectivityType(this->CellToPoint.Connectivity.PrepareForInput(device),
this->CellToPoint.NumIndices.PrepareForInput(device),
this->CellToPoint.IndexOffsets.PrepareForInput(device));
}
void PrintSummary(std::ostream& out) const VTKM_OVERRIDE
VTKM_CONT
void PrintSummary(std::ostream& out) const override
{
out << " CellSetGeneralPermutation of: " << std::endl;
out << "CellSetPermutation of: " << std::endl;
this->FullCellSet.PrintSummary(out);
out << "Permutation Array: " << std::endl;
vtkm::cont::printSummary_ArrayHandle(this->ValidCellIds, out);
}
private:
using CellToPointConnectivity = vtkm::cont::internal::ConnectivityExplicitInternals<
typename ArrayHandleConstant<vtkm::UInt8>::StorageTag>;
PermutationArrayHandleType ValidCellIds;
OriginalCellSet FullCellSet;
OriginalCellSetType FullCellSet;
mutable CellToPointConnectivity CellToPoint;
};
} //namespace internal
#ifndef VTKM_DEFAULT_CELLSET_PERMUTATION_STORAGE_TAG
#define VTKM_DEFAULT_CELLSET_PERMUTATION_STORAGE_TAG VTKM_DEFAULT_STORAGE_TAG
#endif
template <typename OriginalCellSet,
typename PermutationArrayHandleType =
vtkm::cont::ArrayHandle<vtkm::Id, VTKM_DEFAULT_CELLSET_PERMUTATION_STORAGE_TAG>>
class CellSetPermutation
: public vtkm::cont::internal::CellSetGeneralPermutation<OriginalCellSet,
PermutationArrayHandleType>
template <typename CellSetType,
typename PermutationArrayHandleType1,
typename PermutationArrayHandleType2>
class CellSetPermutation<CellSetPermutation<CellSetType, PermutationArrayHandleType1>,
PermutationArrayHandleType2>
: public internal::CellSetPermutationTraits<
CellSetPermutation<CellSetPermutation<CellSetType, PermutationArrayHandleType1>,
PermutationArrayHandleType2>>::Superclass
{
VTKM_IS_CELL_SET(OriginalCellSet);
VTKM_IS_ARRAY_HANDLE(PermutationArrayHandleType);
typedef typename vtkm::cont::internal::CellSetGeneralPermutation<OriginalCellSet,
PermutationArrayHandleType>
ParentType;
private:
using Superclass = typename internal::CellSetPermutationTraits<CellSetPermutation>::Superclass;
public:
VTKM_CONT
CellSetPermutation(const PermutationArrayHandleType& validCellIds,
const OriginalCellSet& cellset,
CellSetPermutation(const PermutationArrayHandleType2& validCellIds,
const CellSetPermutation<CellSetType, PermutationArrayHandleType1>& cellset,
const std::string& name = std::string())
: ParentType(validCellIds, cellset, name)
: Superclass(vtkm::cont::make_ArrayHandlePermutation(validCellIds, cellset.GetValidCellIds()),
cellset.GetFullCellSet(),
name)
{
}
VTKM_CONT
CellSetPermutation(const std::string& name = std::string())
: ParentType(name)
: Superclass(name)
{
}
VTKM_CONT
CellSetPermutation<OriginalCellSet, PermutationArrayHandleType>& operator=(
const CellSetPermutation<OriginalCellSet, PermutationArrayHandleType>& src)
void Fill(const PermutationArrayHandleType2& validCellIds,
const CellSetPermutation<CellSetType, PermutationArrayHandleType1>& cellset)
{
ParentType::operator=(src);
return *this;
this->ValidCellIds = make_ArrayHandlePermutation(validCellIds, cellset.GetValidCellIds());
this->FullCellSet = cellset.GetFullCellSet();
}
};

@ -48,6 +48,7 @@ set(unit_tests
UnitTestArrayHandleConcatenate.cxx
UnitTestArrayPortalToIterators.cxx
UnitTestCellSetExplicit.cxx
UnitTestCellSetPermutation.cxx
UnitTestContTesting.cxx
UnitTestDataSetBuilderExplicit.cxx
UnitTestDataSetBuilderRectilinear.cxx

@ -0,0 +1,195 @@
//============================================================================
// Copyright (c) Kitware, Inc.
// All rights reserved.
// See LICENSE.txt for details.
// This software is distributed WITHOUT ANY WARRANTY; without even
// the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
// PURPOSE. See the above copyright notice for more information.
//
// Copyright 2017 Sandia Corporation.
// Copyright 2017 UT-Battelle, LLC.
// Copyright 2017 Los Alamos National Security.
//
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
//
// Under the terms of Contract DE-AC52-06NA25396 with Los Alamos National
// Laboratory (LANL), the U.S. Government retains certain rights in
// this software.
//============================================================================
#include <vtkm/cont/CellSetPermutation.h>
#include <vtkm/cont/ArrayHandle.h>
#include <vtkm/cont/ArrayHandleCounting.h>
#include <vtkm/cont/DeviceAdapterAlgorithm.h>
#include <vtkm/cont/testing/MakeTestDataSet.h>
#include <vtkm/cont/testing/Testing.h>
#include <vtkm/worklet/DispatcherMapTopology.h>
#include <vtkm/worklet/WorkletMapTopology.h>
namespace
{
struct WorkletPointToCell : public vtkm::worklet::WorkletMapPointToCell
{
typedef void ControlSignature(CellSetIn cellset, FieldOutCell<IdType> numPoints);
typedef void ExecutionSignature(PointIndices, _2);
using InputDomain = _1;
template <typename PointIndicesType>
VTKM_EXEC void operator()(const PointIndicesType& pointIndices, vtkm::Id& numPoints) const
{
numPoints = pointIndices.GetNumberOfComponents();
}
};
struct WorkletCellToPoint : public vtkm::worklet::WorkletMapCellToPoint
{
typedef void ControlSignature(CellSetIn cellset, FieldOutPoint<IdType> numCells);
typedef void ExecutionSignature(CellIndices, _2);
using InputDomain = _1;
template <typename CellIndicesType>
VTKM_EXEC void operator()(const CellIndicesType& cellIndices, vtkm::Id& numCells) const
{
numCells = cellIndices.GetNumberOfComponents();
}
};
struct CellsOfPoint : public vtkm::worklet::WorkletMapCellToPoint
{
typedef void ControlSignature(CellSetIn cellset,
FieldInPoint<IdType> offset,
WholeArrayOut<IdType> cellIds);
typedef void ExecutionSignature(CellIndices, _2, _3);
using InputDomain = _1;
template <typename CellIndicesType, typename CellIdsPortal>
VTKM_EXEC void operator()(const CellIndicesType& cellIndices,
vtkm::Id offset,
const CellIdsPortal& out) const
{
vtkm::IdComponent count = cellIndices.GetNumberOfComponents();
for (vtkm::IdComponent i = 0; i < count; ++i)
{
out.Set(offset++, cellIndices[i]);
}
}
};
template <typename CellSetType, typename PermutationArrayHandleType>
std::vector<vtkm::Id> ComputeCellToPointExpected(const CellSetType& cellset,
const PermutationArrayHandleType& permutation)
{
vtkm::cont::ArrayHandle<vtkm::Id> numIndices;
vtkm::worklet::DispatcherMapTopology<WorkletCellToPoint>().Invoke(cellset, numIndices);
std::cout << "\n";
vtkm::cont::ArrayHandle<vtkm::Id> indexOffsets;
vtkm::Id connectivityLength =
vtkm::cont::DeviceAdapterAlgorithm<VTKM_DEFAULT_DEVICE_ADAPTER_TAG>::ScanExclusive(
numIndices, indexOffsets);
vtkm::cont::ArrayHandle<vtkm::Id> connectivity;
connectivity.Allocate(connectivityLength);
vtkm::worklet::DispatcherMapTopology<CellsOfPoint>().Invoke(cellset, indexOffsets, connectivity);
std::vector<bool> permutationMask(static_cast<std::size_t>(cellset.GetNumberOfCells()), false);
for (vtkm::Id i = 0; i < permutation.GetNumberOfValues(); ++i)
{
permutationMask[static_cast<std::size_t>(permutation.GetPortalConstControl().Get(i))] = true;
}
vtkm::Id numberOfPoints = cellset.GetNumberOfPoints();
std::vector<vtkm::Id> expected(static_cast<std::size_t>(numberOfPoints), 0);
for (vtkm::Id i = 0; i < numberOfPoints; ++i)
{
vtkm::Id offset = indexOffsets.GetPortalConstControl().Get(i);
vtkm::Id count = numIndices.GetPortalConstControl().Get(i);
for (vtkm::Id j = 0; j < count; ++j)
{
vtkm::Id cellId = connectivity.GetPortalConstControl().Get(offset++);
if (permutationMask[static_cast<std::size_t>(cellId)])
{
++expected[static_cast<std::size_t>(i)];
}
}
}
return expected;
}
template <typename CellSetType>
vtkm::cont::CellSetPermutation<CellSetType, vtkm::cont::ArrayHandleCounting<vtkm::Id>> TestCellSet(
const CellSetType& cellset)
{
vtkm::Id numberOfCells = cellset.GetNumberOfCells() / 2;
vtkm::cont::ArrayHandleCounting<vtkm::Id> permutation(0, 2, numberOfCells);
auto cs = vtkm::cont::make_CellSetPermutation(permutation, cellset);
vtkm::cont::ArrayHandle<vtkm::Id> result;
std::cout << "\t\tTesting PointToCell\n";
vtkm::worklet::DispatcherMapTopology<WorkletPointToCell>().Invoke(cs, result);
VTKM_TEST_ASSERT(result.GetNumberOfValues() == numberOfCells,
"result length not equal to number of cells");
for (vtkm::Id i = 0; i < result.GetNumberOfValues(); ++i)
{
VTKM_TEST_ASSERT(result.GetPortalConstControl().Get(i) ==
cellset.GetNumberOfPointsInCell(permutation.GetPortalConstControl().Get(i)),
"incorrect result");
}
std::cout << "\t\tTesting CellToPoint\n";
vtkm::worklet::DispatcherMapTopology<WorkletCellToPoint>().Invoke(cs, result);
VTKM_TEST_ASSERT(result.GetNumberOfValues() == cellset.GetNumberOfPoints(),
"result length not equal to number of points");
auto expected = ComputeCellToPointExpected(cellset, permutation);
for (vtkm::Id i = 0; i < result.GetNumberOfValues(); ++i)
{
VTKM_TEST_ASSERT(result.GetPortalConstControl().Get(i) == expected[static_cast<std::size_t>(i)],
"incorrect result");
}
return cs;
}
template <typename CellSetType>
void RunTests(const CellSetType& cellset)
{
std::cout << "\tTesting CellSetPermutation:\n";
auto p1 = TestCellSet(cellset);
std::cout << "\tTesting CellSetPermutation of CellSetPermutation:\n";
TestCellSet(p1);
std::cout << "----------------------------------------------------------\n";
}
void TestCellSetPermutation()
{
vtkm::cont::DataSet dataset;
vtkm::cont::testing::MakeTestDataSet maker;
std::cout << "Testing CellSetStructured<2>\n";
dataset = maker.Make2DUniformDataSet1();
RunTests(dataset.GetCellSet().Cast<vtkm::cont::CellSetStructured<2>>());
std::cout << "Testing CellSetStructured<3>\n";
dataset = maker.Make3DUniformDataSet1();
RunTests(dataset.GetCellSet().Cast<vtkm::cont::CellSetStructured<3>>());
std::cout << "Testing CellSetExplicit\n";
dataset = maker.Make3DExplicitDataSetPolygonal();
RunTests(dataset.GetCellSet().Cast<vtkm::cont::CellSetExplicit<>>());
std::cout << "Testing CellSetSingleType\n";
dataset = maker.Make3DExplicitDataSetCowNose();
RunTests(dataset.GetCellSet().Cast<vtkm::cont::CellSetSingleType<>>());
}
} // anonymous namespace
int UnitTestCellSetPermutation(int, char* [])
{
return vtkm::cont::testing::Testing::Run(TestCellSetPermutation);
}

@ -21,9 +21,10 @@
#ifndef vtk_m_exec_ConnectivityPermuted_h
#define vtk_m_exec_ConnectivityPermuted_h
#include <vtkm/CellShape.h>
#include <vtkm/TopologyElementTag.h>
#include <vtkm/Types.h>
#include <vtkm/exec/ConnectivityStructured.h>
#include <vtkm/VecFromPortal.h>
namespace vtkm
{
@ -31,28 +32,28 @@ namespace exec
{
template <typename PermutationPortal, typename OriginalConnectivity>
class ConnectivityPermuted
class ConnectivityPermutedPointToCell
{
public:
typedef typename OriginalConnectivity::SchedulingRangeType SchedulingRangeType;
VTKM_SUPPRESS_EXEC_WARNINGS
VTKM_EXEC_CONT
ConnectivityPermuted()
ConnectivityPermutedPointToCell()
: Portal()
, Connectivity()
{
}
VTKM_EXEC_CONT
ConnectivityPermuted(const PermutationPortal& portal, const OriginalConnectivity& src)
ConnectivityPermutedPointToCell(const PermutationPortal& portal, const OriginalConnectivity& src)
: Portal(portal)
, Connectivity(src)
{
}
VTKM_EXEC_CONT
ConnectivityPermuted(const ConnectivityPermuted& src)
ConnectivityPermutedPointToCell(const ConnectivityPermutedPointToCell& src)
: Portal(src.Portal)
, Connectivity(src.Connectivity)
{
@ -87,6 +88,47 @@ public:
PermutationPortal Portal;
OriginalConnectivity Connectivity;
};
template <typename ConnectivityPortalType,
typename NumIndicesPortalType,
typename IndexOffsetPortalType>
class ConnectivityPermutedCellToPoint
{
public:
using SchedulingRangeType = vtkm::Id;
using IndicesType = vtkm::VecFromPortal<ConnectivityPortalType>;
using CellShapeTag = vtkm::CellShapeTagVertex;
ConnectivityPermutedCellToPoint() = default;
ConnectivityPermutedCellToPoint(const ConnectivityPortalType& connectivity,
const NumIndicesPortalType& numIndices,
const IndexOffsetPortalType& indexOffset)
: Connectivity(connectivity)
, NumIndices(numIndices)
, IndexOffset(indexOffset)
{
}
VTKM_EXEC
SchedulingRangeType GetNumberOfElements() const { return this->NumIndices.GetNumberOfValues(); }
VTKM_EXEC CellShapeTag GetCellShape(vtkm::Id) const { return CellShapeTag(); }
VTKM_EXEC
vtkm::IdComponent GetNumberOfIndices(vtkm::Id index) const { return this->NumIndices.Get(index); }
VTKM_EXEC IndicesType GetIndices(vtkm::Id index) const
{
return IndicesType(
this->Connectivity, this->NumIndices.Get(index), this->IndexOffset.Get(index));
}
private:
ConnectivityPortalType Connectivity;
NumIndicesPortalType NumIndices;
IndexOffsetPortalType IndexOffset;
};
}
} // namespace vtkm::exec

@ -159,14 +159,15 @@ struct FetchArrayTopologyMapInImplementation<
template <typename PermutationPortal, vtkm::IdComponent NumDimensions>
struct FetchArrayTopologyMapInImplementation<
vtkm::exec::ConnectivityPermuted<PermutationPortal,
vtkm::exec::ConnectivityStructured<vtkm::TopologyElementTagPoint,
vtkm::TopologyElementTagCell,
NumDimensions>>,
vtkm::exec::ConnectivityPermutedPointToCell<
PermutationPortal,
vtkm::exec::ConnectivityStructured<vtkm::TopologyElementTagPoint,
vtkm::TopologyElementTagCell,
NumDimensions>>,
vtkm::internal::ArrayPortalUniformPointCoordinates>
{
typedef vtkm::exec::ConnectivityPermuted<
typedef vtkm::exec::ConnectivityPermutedPointToCell<
PermutationPortal,
vtkm::exec::ConnectivityStructured<vtkm::TopologyElementTagPoint,
vtkm::TopologyElementTagCell,

@ -347,18 +347,20 @@ private:
};
// Specialization for permuted structured connectivity types.
template <typename PermutationPortal,
typename FromTopology,
typename ToTopology,
vtkm::IdComponent Dimension>
class ThreadIndicesTopologyMap<vtkm::exec::ConnectivityPermuted<
template <typename PermutationPortal, vtkm::IdComponent Dimension>
class ThreadIndicesTopologyMap<vtkm::exec::ConnectivityPermutedPointToCell<
PermutationPortal,
vtkm::exec::ConnectivityStructured<FromTopology, ToTopology, Dimension>>>
vtkm::exec::
ConnectivityStructured<vtkm::TopologyElementTagPoint, vtkm::TopologyElementTagCell, Dimension>>>
{
using PermutedConnectivityType = vtkm::exec::ConnectivityPermuted<
using PermutedConnectivityType = vtkm::exec::ConnectivityPermutedPointToCell<
PermutationPortal,
vtkm::exec::ConnectivityStructured<FromTopology, ToTopology, Dimension>>;
using ConnectivityType = vtkm::exec::ConnectivityStructured<FromTopology, ToTopology, Dimension>;
vtkm::exec::ConnectivityStructured<vtkm::TopologyElementTagPoint,
vtkm::TopologyElementTagCell,
Dimension>>;
using ConnectivityType = vtkm::exec::ConnectivityStructured<vtkm::TopologyElementTagPoint,
vtkm::TopologyElementTagCell,
Dimension>;
public:
using IndicesFromType = typename ConnectivityType::IndicesType;