From 5a22c135a9cdbfa07380562ceb920d84c106983c Mon Sep 17 00:00:00 2001 From: Kenneth Moreland <morelandkd@ornl.gov> Date: Fri, 15 Dec 2023 10:13:09 -0700 Subject: [PATCH] Fix issues with VTK to VTK-m array conversion Change the VTK to VTK-m array conversion routines to use `ArrayHandleRuntimeVec` and `ArrayHandleRecombineVec`. These are new features of VTK-m that allow you to specify an array with the tuple size specified at runtime. This change improves three specific things. * Fixes a bug when importing an array of "odd" tuple size (not 1, 2, 3, 4, 6, or 9). It was creating arrays of size one less than the actual size. * Avoids using `ArrayHandleGroupVecVariable`, which is supported by fewer VTK-m filters. * The VTK-m `ArrayHandle` now manages a reference back to the VTK array, so the `ArrayHandle` will continue to work even if the original VTK array is "deleted." This makes the code safer. * Unifies the implementation of the array conversion among number of components to avoid issues with surprise tuple sizes. --- .../Vtkm/Core/vtkmlib/DataArrayConverters.h | 98 +++++++++++++++++++ .../Vtkm/Core/vtkmlib/DataArrayConverters.hxx | 30 ------ .../DataModel/vtkmlib/CellSetConverters.cxx | 15 +-- .../DataModel/vtkmlib/DataSetConverters.cxx | 8 +- Accelerators/Vtkm/Filters/vtkmGradient.cxx | 3 +- Accelerators/Vtkm/Filters/vtkmThreshold.cxx | 3 +- .../dev/corrected-vtk-to-vtkm-arrays.md | 17 ++++ ThirdParty/vtkm/vtkvtkm/vtk-m | 2 +- 8 files changed, 126 insertions(+), 50 deletions(-) create mode 100644 Documentation/release/dev/corrected-vtk-to-vtkm-arrays.md diff --git a/Accelerators/Vtkm/Core/vtkmlib/DataArrayConverters.h b/Accelerators/Vtkm/Core/vtkmlib/DataArrayConverters.h index b6b2ce80626..6e8ec1df50b 100644 --- a/Accelerators/Vtkm/Core/vtkmlib/DataArrayConverters.h +++ b/Accelerators/Vtkm/Core/vtkmlib/DataArrayConverters.h @@ -12,11 +12,19 @@ #include "vtkAOSDataArrayTemplate.h" #include "vtkSOADataArrayTemplate.h" +#include "vtkLogger.h" + +#include <vtkm/cont/ArrayExtractComponent.h> +#include <vtkm/cont/ArrayHandleBasic.h> +#include <vtkm/cont/ArrayHandleRecombineVec.h> +#include <vtkm/cont/ArrayHandleRuntimeVec.h> #include <vtkm/cont/ArrayHandleSOA.h> +#include <vtkm/cont/ArrayHandleStride.h> #include <vtkm/cont/Field.h> #include <vtkm/cont/UnknownArrayHandle.h> #include <type_traits> // for std::underlying_type +#include <utility> // for std::pair namespace vtkm { @@ -44,6 +52,93 @@ inline static const char* NoNameVTKFieldName() return name; } +template <typename T> +vtkm::cont::ArrayHandleBasic<T> vtkAOSDataArrayToFlatArrayHandle(vtkAOSDataArrayTemplate<T>* input) +{ + // Register a reference to the input here to make sure the array cannot + // be deleted before the `ArrayHandle` is done with it. (Note that you + // will still get problems if the `vtkAOSDataArrayTemplate` gets resized. + input->Register(nullptr); + + auto deleter = [](void* container) { + vtkAOSDataArrayTemplate<T>* vtkArray = reinterpret_cast<vtkAOSDataArrayTemplate<T>*>(container); + vtkArray->UnRegister(nullptr); + }; + auto reallocator = [](void*& memory, void*& container, vtkm::BufferSizeType oldSize, + vtkm::BufferSizeType newSize) { + vtkAOSDataArrayTemplate<T>* vtkArray = reinterpret_cast<vtkAOSDataArrayTemplate<T>*>(container); + if ((vtkArray->GetVoidPointer(0) != memory) || (vtkArray->GetNumberOfValues() != oldSize)) + { + vtkLog(ERROR, + "Dangerous inconsistency found between pointers for VTK and VTK-m. " + "Was the VTK array resized outside of VTK-m?"); + } + vtkArray->SetNumberOfValues(newSize); + memory = vtkArray->GetVoidPointer(0); + }; + + return vtkm::cont::ArrayHandleBasic<T>( + input->GetPointer(0), input, input->GetNumberOfValues(), deleter, reallocator); +} + +template <typename T> +vtkm::cont::ArrayHandleBasic<T> vtkSOADataArrayToComponentArrayHandle( + vtkSOADataArrayTemplate<T>* input, int componentIndex) +{ + // Register for each component (as each will have the deleter call to + // unregister). + input->Register(nullptr); + + using ContainerPair = std::pair<vtkSOADataArrayTemplate<T>*, int>; + ContainerPair* componentInput = new ContainerPair(input, componentIndex); + + auto deleter = [](void* container) { + ContainerPair* containerPair = reinterpret_cast<ContainerPair*>(container); + containerPair->first->UnRegister(nullptr); + delete containerPair; + }; + auto reallocator = [](void*& memory, void*& container, vtkm::BufferSizeType vtkNotUsed(oldSize), + vtkm::BufferSizeType newSize) { + ContainerPair* containerPair = reinterpret_cast<ContainerPair*>(container); + containerPair->first->SetNumberOfTuples(newSize); + memory = containerPair->first->GetComponentArrayPointer(containerPair->second); + }; + + return vtkm::cont::ArrayHandleBasic<T>(input->GetComponentArrayPointer(componentIndex), + componentInput, input->GetNumberOfTuples(), deleter, reallocator); +} + +template <typename T> +vtkm::cont::ArrayHandleRuntimeVec<T> vtkDataArrayToArrayHandle(vtkAOSDataArrayTemplate<T>* input) +{ + auto flatArray = vtkAOSDataArrayToFlatArrayHandle(input); + return vtkm::cont::make_ArrayHandleRuntimeVec(input->GetNumberOfComponents(), flatArray); +} + +template <typename T> +vtkm::cont::ArrayHandleRecombineVec<T> vtkDataArrayToArrayHandle(vtkSOADataArrayTemplate<T>* input) +{ + // Wrap each component array in a basic array handle, convert that to a + // strided array, and then add that as a component to the returned + // recombined vec. + vtkm::cont::ArrayHandleRecombineVec<T> output; + + for (int componentIndex = 0; componentIndex < input->GetNumberOfComponents(); ++componentIndex) + { + auto componentArray = vtkSOADataArrayToComponentArrayHandle(input, componentIndex); + output.AppendComponentArray( + vtkm::cont::ArrayExtractComponent(componentArray, 0, vtkm::CopyFlag::Off)); + } + + return output; +} + +template <typename DataArrayType> +vtkm::cont::UnknownArrayHandle vtkDataArrayToUnknownArrayHandle(DataArrayType* input) +{ + return vtkDataArrayToArrayHandle(input); +} + template <typename DataArrayType, vtkm::IdComponent NumComponents> struct DataArrayToArrayHandle; @@ -55,6 +150,7 @@ struct DataArrayToArrayHandle<vtkAOSDataArrayTemplate<T>, NumComponents> using StorageType = vtkm::cont::internal::Storage<ValueType, vtkm::cont::StorageTagBasic>; using ArrayHandleType = vtkm::cont::ArrayHandle<ValueType, vtkm::cont::StorageTagBasic>; + VTK_DEPRECATED_IN_9_3_0("Use vtkDataArrayToArrayHandle or vtkAOSDataArrayToFlatArrayHandle.") static ArrayHandleType Wrap(vtkAOSDataArrayTemplate<T>* input) { return vtkm::cont::make_ArrayHandle(reinterpret_cast<ValueType*>(input->GetPointer(0)), @@ -69,6 +165,7 @@ struct DataArrayToArrayHandle<vtkSOADataArrayTemplate<T>, NumComponents> using StorageType = vtkm::cont::internal::Storage<ValueType, vtkm::cont::StorageTagSOA>; using ArrayHandleType = vtkm::cont::ArrayHandle<ValueType, vtkm::cont::StorageTagSOA>; + VTK_DEPRECATED_IN_9_3_0("Use vtkDataArrayToArrayHandle or vtkSOADataArrayToComponentArrayHandle.") static ArrayHandleType Wrap(vtkSOADataArrayTemplate<T>* input) { vtkm::Id numValues = input->GetNumberOfTuples(); @@ -90,6 +187,7 @@ struct DataArrayToArrayHandle<vtkSOADataArrayTemplate<T>, 1> using StorageType = vtkm::cont::internal::Storage<T, vtkm::cont::StorageTagBasic>; using ArrayHandleType = vtkm::cont::ArrayHandle<T, vtkm::cont::StorageTagBasic>; + VTK_DEPRECATED_IN_9_3_0("Use vtkDataArrayToArrayHandle or vtkSOADataArrayToComponentArrayHandle.") static ArrayHandleType Wrap(vtkSOADataArrayTemplate<T>* input) { return vtkm::cont::make_ArrayHandle( diff --git a/Accelerators/Vtkm/Core/vtkmlib/DataArrayConverters.hxx b/Accelerators/Vtkm/Core/vtkmlib/DataArrayConverters.hxx index a730d162b27..9805822cd4b 100644 --- a/Accelerators/Vtkm/Core/vtkmlib/DataArrayConverters.hxx +++ b/Accelerators/Vtkm/Core/vtkmlib/DataArrayConverters.hxx @@ -16,36 +16,6 @@ namespace tovtkm { VTK_ABI_NAMESPACE_BEGIN -template <typename DataArrayType> -vtkm::cont::UnknownArrayHandle vtkDataArrayToUnknownArrayHandle(DataArrayType* input) -{ - int numComps = input->GetNumberOfComponents(); - switch (numComps) - { - case 1: - return vtkm::cont::UnknownArrayHandle(DataArrayToArrayHandle<DataArrayType, 1>::Wrap(input)); - case 2: - return vtkm::cont::UnknownArrayHandle(DataArrayToArrayHandle<DataArrayType, 2>::Wrap(input)); - case 3: - return vtkm::cont::UnknownArrayHandle(DataArrayToArrayHandle<DataArrayType, 3>::Wrap(input)); - case 4: - return vtkm::cont::UnknownArrayHandle(DataArrayToArrayHandle<DataArrayType, 4>::Wrap(input)); - case 6: - return vtkm::cont::UnknownArrayHandle(DataArrayToArrayHandle<DataArrayType, 6>::Wrap(input)); - case 9: - return vtkm::cont::UnknownArrayHandle(DataArrayToArrayHandle<DataArrayType, 9>::Wrap(input)); - default: - { - vtkm::Id numTuples = input->GetNumberOfTuples(); - auto subHandle = DataArrayToArrayHandle<DataArrayType, 1>::Wrap(input); - auto offsets = - vtkm::cont::ArrayHandleCounting<vtkm::Id>(vtkm::Id(0), vtkm::Id(numComps), numTuples); - auto handle = vtkm::cont::make_ArrayHandleGroupVecVariable(subHandle, offsets); - return vtkm::cont::UnknownArrayHandle(handle); - } - } -} - template <typename DataArrayType> vtkm::cont::Field ConvertPointField(DataArrayType* input) { diff --git a/Accelerators/Vtkm/DataModel/vtkmlib/CellSetConverters.cxx b/Accelerators/Vtkm/DataModel/vtkmlib/CellSetConverters.cxx index 112627c9d89..2b7001e8f13 100644 --- a/Accelerators/Vtkm/DataModel/vtkmlib/CellSetConverters.cxx +++ b/Accelerators/Vtkm/DataModel/vtkmlib/CellSetConverters.cxx @@ -55,12 +55,10 @@ struct BuildSingleTypeCellSetVisitor CellStateT& state, vtkm::UInt8 cellType, vtkm::IdComponent cellSize, vtkIdType numPoints) { using VTKIdT = typename CellStateT::ValueType; // might not be vtkIdType... - using VTKArrayT = vtkAOSDataArrayTemplate<VTKIdT>; static constexpr bool IsVtkmIdType = std::is_same<VTKIdT, vtkm::Id>::value; // Construct an arrayhandle that holds the connectivity array - using DirectConverter = tovtkm::DataArrayToArrayHandle<VTKArrayT, 1>; - auto connHandleDirect = DirectConverter::Wrap(state.GetConnectivity()); + auto connHandleDirect = tovtkm::vtkAOSDataArrayToFlatArrayHandle(state.GetConnectivity()); // Cast if necessary: auto connHandle = IsVtkmIdType ? connHandleDirect @@ -164,13 +162,11 @@ struct BuildExplicitCellSetVisitor const vtkm::cont::ArrayHandle<vtkm::UInt8, S>& shapes, vtkm::Id numPoints) const { using VTKIdT = typename CellStateT::ValueType; // might not be vtkIdType... - using VTKArrayT = vtkAOSDataArrayTemplate<VTKIdT>; static constexpr bool IsVtkmIdType = std::is_same<VTKIdT, vtkm::Id>::value; // Construct arrayhandles to hold the arrays - using DirectConverter = tovtkm::DataArrayToArrayHandle<VTKArrayT, 1>; - auto offsetsHandleDirect = DirectConverter::Wrap(state.GetOffsets()); - auto connHandleDirect = DirectConverter::Wrap(state.GetConnectivity()); + auto offsetsHandleDirect = tovtkm::vtkAOSDataArrayToFlatArrayHandle(state.GetOffsets()); + auto connHandleDirect = tovtkm::vtkAOSDataArrayToFlatArrayHandle(state.GetConnectivity()); // Cast if necessary: auto connHandle = IsVtkmIdType ? connHandleDirect @@ -207,10 +203,7 @@ struct SupportedCellShape vtkm::cont::UnknownCellSet Convert( vtkUnsignedCharArray* types, vtkCellArray* cells, vtkIdType numberOfPoints) { - using ShapeArrayType = vtkAOSDataArrayTemplate<vtkm::UInt8>; - using ShapeConverter = tovtkm::DataArrayToArrayHandle<ShapeArrayType, 1>; - - auto shapes = ShapeConverter::Wrap(types); + auto shapes = tovtkm::vtkAOSDataArrayToFlatArrayHandle(types); if (!vtkm::cont::Algorithm::Reduce( vtkm::cont::make_ArrayHandleTransform(shapes, SupportedCellShape{}), true, vtkm::LogicalAnd())) diff --git a/Accelerators/Vtkm/DataModel/vtkmlib/DataSetConverters.cxx b/Accelerators/Vtkm/DataModel/vtkmlib/DataSetConverters.cxx index ee3bb566829..f3fcc84ceab 100644 --- a/Accelerators/Vtkm/DataModel/vtkmlib/DataSetConverters.cxx +++ b/Accelerators/Vtkm/DataModel/vtkmlib/DataSetConverters.cxx @@ -49,7 +49,7 @@ vtkm::cont::CoordinateSystem deduce_container(vtkPoints* points) vtkAOSDataArrayTemplate<T>* typedIn = vtkAOSDataArrayTemplate<T>::FastDownCast(points->GetData()); if (typedIn) { - auto p = DataArrayToArrayHandle<vtkAOSDataArrayTemplate<T>, 3>::Wrap(typedIn); + auto p = vtkDataArrayToArrayHandle(typedIn); return vtkm::cont::CoordinateSystem("coords", p); } @@ -57,7 +57,7 @@ vtkm::cont::CoordinateSystem deduce_container(vtkPoints* points) vtkSOADataArrayTemplate<T>::FastDownCast(points->GetData()); if (typedIn2) { - auto p = DataArrayToArrayHandle<vtkSOADataArrayTemplate<T>, 3>::Wrap(typedIn2); + auto p = vtkDataArrayToArrayHandle(typedIn2); return vtkm::cont::CoordinateSystem("coords", p); } @@ -147,14 +147,14 @@ vtkm::cont::CoordinateSystem ConvertRectilinearPoints( vtkAOSDataArrayTemplate<T>* aos = vtkAOSDataArrayTemplate<T>::FastDownCast(vtkCompArrays[i]); if (aos) { - vtkmCompArrays[i] = DataArrayToArrayHandle<vtkAOSDataArrayTemplate<T>, 1>::Wrap(aos); + vtkmCompArrays[i] = vtkAOSDataArrayToFlatArrayHandle(aos); continue; } vtkSOADataArrayTemplate<T>* soa = vtkSOADataArrayTemplate<T>::FastDownCast(vtkCompArrays[i]); if (soa) { - vtkmCompArrays[i] = DataArrayToArrayHandle<vtkSOADataArrayTemplate<T>, 1>::Wrap(soa); + vtkmCompArrays[i] = vtkSOADataArrayToComponentArrayHandle(soa, 0); continue; } diff --git a/Accelerators/Vtkm/Filters/vtkmGradient.cxx b/Accelerators/Vtkm/Filters/vtkmGradient.cxx index 73e69e56262..d2efdf0f581 100644 --- a/Accelerators/Vtkm/Filters/vtkmGradient.cxx +++ b/Accelerators/Vtkm/Filters/vtkmGradient.cxx @@ -59,8 +59,7 @@ private: inline bool HasGhostFlagsSet(vtkUnsignedCharArray* ghostArray, int flags) { - auto ah = - tovtkm::DataArrayToArrayHandle<vtkAOSDataArrayTemplate<unsigned char>, 1>::Wrap(ghostArray); + auto ah = tovtkm::vtkAOSDataArrayToFlatArrayHandle(ghostArray); int result = vtkm::cont::Algorithm::Reduce( vtkm::cont::make_ArrayHandleTransform(ah, MaskBits(flags)), 0, vtkm::LogicalOr()); return result; diff --git a/Accelerators/Vtkm/Filters/vtkmThreshold.cxx b/Accelerators/Vtkm/Filters/vtkmThreshold.cxx index fad613e6939..831623f4ab3 100644 --- a/Accelerators/Vtkm/Filters/vtkmThreshold.cxx +++ b/Accelerators/Vtkm/Filters/vtkmThreshold.cxx @@ -56,8 +56,7 @@ inline bool HasGhostFlagsSet(vtkUnsignedCharArray* ghostArray, int flags) return false; } - auto ah = - tovtkm::DataArrayToArrayHandle<vtkAOSDataArrayTemplate<unsigned char>, 1>::Wrap(ghostArray); + auto ah = tovtkm::vtkAOSDataArrayToFlatArrayHandle(ghostArray); int result = vtkm::cont::Algorithm::Reduce( vtkm::cont::make_ArrayHandleTransform(ah, MaskBits(flags)), 0, vtkm::LogicalOr()); return result; diff --git a/Documentation/release/dev/corrected-vtk-to-vtkm-arrays.md b/Documentation/release/dev/corrected-vtk-to-vtkm-arrays.md new file mode 100644 index 00000000000..d124bda3d74 --- /dev/null +++ b/Documentation/release/dev/corrected-vtk-to-vtkm-arrays.md @@ -0,0 +1,17 @@ +Fix issues with VTK to VTK-m array conversion + +Change the VTK to VTK-m array conversion routines to use +`ArrayHandleRuntimeVec` and `ArrayHandleRecombineVec`. These are new +features of VTK-m that allow you to specify an array with the tuple size +specified at runtime. This change improves three specific things. + +* Fixes a bug when importing an array of "odd" tuple size (not 1, 2, + 3, 4, 6, or 9). It was creating arrays of size one less than the + actual size. +* Avoids using `ArrayHandleGroupVecVariable`, which is supported by + fewer VTK-m filters. +* The VTK-m `ArrayHandle` now manages a reference back to the VTK + array, so the `ArrayHandle` will continue to work even if the + original VTK array is "deleted." This makes the code safer. +* Unifies the implementation of the array conversion among number of + components to avoid issues with surprise tuple sizes. diff --git a/ThirdParty/vtkm/vtkvtkm/vtk-m b/ThirdParty/vtkm/vtkvtkm/vtk-m index a057f62e756..3c9249871a5 160000 --- a/ThirdParty/vtkm/vtkvtkm/vtk-m +++ b/ThirdParty/vtkm/vtkvtkm/vtk-m @@ -1 +1 @@ -Subproject commit a057f62e756efc43095e72c5813aaaf2dea36ebb +Subproject commit 3c9249871a5458cac01bbfe83307232aadd6bed5 -- GitLab