//============================================================================= // // 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 2012 Sandia Corporation. // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, // the U.S. Government retains certain rights in this software. // //============================================================================= #include "CellSetConverters.h" #include "ArrayConverters.hxx" #include "DataSetConverters.h" #include "vtkmFilterPolicy.h" #include #include #include #include #include #include #include #include #include #include "vtkCellArray.h" #include "vtkCellType.h" #include "vtkIdTypeArray.h" #include "vtkNew.h" #include "vtkUnsignedCharArray.h" #if defined(VTKM_CUDA) #define FUNC_SCOPE __device__ #else #define FUNC_SCOPE #endif namespace tovtkm { namespace { struct ReorderHex : vtkm::worklet::WorkletMapField { using ControlSignature = void(FieldInOut); FUNC_SCOPE void operator()(vtkm::Vec& indices) const { auto doSwap = [&](vtkm::IdComponent id1, vtkm::IdComponent id2) { const auto t = indices[id1]; indices[id1] = indices[id2]; indices[id2] = t; }; doSwap(2, 3); doSwap(6, 7); } }; struct BuildSingleTypeCellSetVisitor { template vtkm::cont::UnknownCellSet operator()( CellStateT& state, vtkm::UInt8 cellType, vtkm::IdComponent cellSize, vtkIdType numPoints) { using VTKIdT = typename CellStateT::ValueType; // might not be vtkIdType... using VTKArrayT = vtkAOSDataArrayTemplate; static constexpr bool IsVtkmIdType = std::is_same::value; // Construct an arrayhandle that holds the connectivity array using DirectConverter = tovtkm::DataArrayToArrayHandle; auto connHandleDirect = DirectConverter::Wrap(state.GetConnectivity()); // Cast if necessary: auto connHandle = IsVtkmIdType ? connHandleDirect : vtkm::cont::make_ArrayHandleCast(connHandleDirect); using ConnHandleType = typename std::decay::type; using ConnStorageTag = typename ConnHandleType::StorageTag; using CellSetType = vtkm::cont::CellSetSingleType; CellSetType cellSet; cellSet.Fill(static_cast(numPoints), cellType, cellSize, connHandle); return cellSet; } }; struct BuildSingleTypeVoxelCellSetVisitor { template vtkm::cont::UnknownCellSet operator()(CellStateT& state, vtkIdType numPoints) { vtkm::cont::ArrayHandle connHandle; { auto* conn = state.GetConnectivity(); const auto* origData = conn->GetPointer(0); const vtkm::Id numIds = conn->GetNumberOfValues(); vtkm::cont::ArrayCopy( vtkm::cont::make_ArrayHandle(origData, numIds, vtkm::CopyFlag::Off), connHandle); // reorder cells from voxel->hex: which only can run on // devices that have shared memory / vtable with the CPU vtkm::cont::ScopedRuntimeDeviceTracker tracker( vtkm::cont::DeviceAdapterTagAny{}, vtkm::cont::RuntimeDeviceTrackerMode::Disable); tracker.ResetDevice(vtkm::cont::DeviceAdapterTagTBB{}); tracker.ResetDevice(vtkm::cont::DeviceAdapterTagOpenMP{}); tracker.ResetDevice(vtkm::cont::DeviceAdapterTagSerial{}); vtkm::cont::Invoker invoke; invoke(ReorderHex{}, vtkm::cont::make_ArrayHandleGroupVec<8>(connHandle)); } using CellSetType = vtkm::cont::CellSetSingleType<>; CellSetType cellSet; cellSet.Fill(numPoints, vtkm::CELL_SHAPE_HEXAHEDRON, 8, connHandle); return cellSet; } }; } // end anon namespace // convert a cell array of a single type to a vtkm CellSetSingleType vtkm::cont::UnknownCellSet ConvertSingleType( vtkCellArray* cells, int cellType, vtkIdType numberOfPoints) { switch (cellType) { case VTK_LINE: return cells->Visit( BuildSingleTypeCellSetVisitor{}, vtkm::CELL_SHAPE_LINE, 2, numberOfPoints); case VTK_HEXAHEDRON: return cells->Visit( BuildSingleTypeCellSetVisitor{}, vtkm::CELL_SHAPE_HEXAHEDRON, 8, numberOfPoints); case VTK_VOXEL: // Note that this is a special case that reorders ids voxel to hex: return cells->Visit(BuildSingleTypeVoxelCellSetVisitor{}, numberOfPoints); case VTK_QUAD: return cells->Visit( BuildSingleTypeCellSetVisitor{}, vtkm::CELL_SHAPE_QUAD, 4, numberOfPoints); case VTK_TETRA: return cells->Visit( BuildSingleTypeCellSetVisitor{}, vtkm::CELL_SHAPE_TETRA, 4, numberOfPoints); case VTK_TRIANGLE: return cells->Visit( BuildSingleTypeCellSetVisitor{}, vtkm::CELL_SHAPE_TRIANGLE, 3, numberOfPoints); case VTK_VERTEX: return cells->Visit( BuildSingleTypeCellSetVisitor{}, vtkm::CELL_SHAPE_VERTEX, 1, numberOfPoints); case VTK_WEDGE: return cells->Visit( BuildSingleTypeCellSetVisitor{}, vtkm::CELL_SHAPE_WEDGE, 6, numberOfPoints); case VTK_PYRAMID: return cells->Visit( BuildSingleTypeCellSetVisitor{}, vtkm::CELL_SHAPE_PYRAMID, 5, numberOfPoints); default: break; } throw vtkm::cont::ErrorBadType("Unsupported VTK cell type in " "CellSetSingleType converter."); } namespace { struct BuildExplicitCellSetVisitor { template vtkm::cont::UnknownCellSet operator()(CellStateT& state, const vtkm::cont::ArrayHandle& shapes, vtkm::Id numPoints) const { using VTKIdT = typename CellStateT::ValueType; // might not be vtkIdType... using VTKArrayT = vtkAOSDataArrayTemplate; static constexpr bool IsVtkmIdType = std::is_same::value; // Construct arrayhandles to hold the arrays using DirectConverter = tovtkm::DataArrayToArrayHandle; auto offsetsHandleDirect = DirectConverter::Wrap(state.GetOffsets()); auto connHandleDirect = DirectConverter::Wrap(state.GetConnectivity()); // Cast if necessary: auto connHandle = IsVtkmIdType ? connHandleDirect : vtkm::cont::make_ArrayHandleCast(connHandleDirect); auto offsetsHandle = IsVtkmIdType ? offsetsHandleDirect : vtkm::cont::make_ArrayHandleCast(offsetsHandleDirect); using ShapesStorageTag = typename std::decay::type::StorageTag; using ConnStorageTag = typename decltype(connHandle)::StorageTag; using OffsetsStorageTag = typename decltype(offsetsHandle)::StorageTag; using CellSetType = vtkm::cont::CellSetExplicit; CellSetType cellSet; cellSet.Fill(numPoints, shapes, connHandle, offsetsHandle); return cellSet; } }; } // end anon namespace // convert a cell array of mixed types to a vtkm CellSetExplicit vtkm::cont::UnknownCellSet Convert( vtkUnsignedCharArray* types, vtkCellArray* cells, vtkIdType numberOfPoints) { using ShapeArrayType = vtkAOSDataArrayTemplate; using ShapeConverter = tovtkm::DataArrayToArrayHandle; return cells->Visit(BuildExplicitCellSetVisitor{}, ShapeConverter::Wrap(types), numberOfPoints); } } // namespace tovtkm namespace fromvtkm { bool Convert(const vtkm::cont::UnknownCellSet& toConvert, vtkCellArray* cells, vtkUnsignedCharArray* typesArray) { const auto* cellset = toConvert.GetCellSetBase(); // small hack as we can't compute properly the number of cells // instead we will pre-allocate and than shrink const vtkm::Id numCells = cellset->GetNumberOfCells(); const vtkm::Id maxSize = numCells * 8; // largest cell type is hex // TODO this could steal the guts out of explicit cellsets as a future // no-copy optimization. vtkNew offsetsArray; offsetsArray->SetNumberOfTuples(static_cast(numCells + 1)); vtkNew connArray; connArray->SetNumberOfTuples(static_cast(maxSize)); if (typesArray) { typesArray->SetNumberOfComponents(1); typesArray->SetNumberOfTuples(static_cast(numCells)); } vtkIdType* connIter = connArray->GetPointer(0); const vtkIdType* connBegin = connIter; for (vtkm::Id cellId = 0; cellId < numCells; ++cellId) { const vtkIdType vtkCellId = static_cast(cellId); const vtkm::Id npts = cellset->GetNumberOfPointsInCell(cellId); assert(npts <= 8 && "Initial allocation assumes no more than 8 pts/cell."); const vtkIdType offset = static_cast(connIter - connBegin); offsetsArray->SetValue(vtkCellId, offset); cellset->GetCellPointIds(cellId, connIter); connIter += npts; if (typesArray) { typesArray->SetValue(vtkCellId, cellset->GetCellShape(cellId)); } } const vtkIdType connSize = static_cast(connIter - connBegin); offsetsArray->SetValue(static_cast(numCells), connSize); connArray->Resize(connSize); cells->SetData(offsetsArray, connArray); return true; } } // namespace fromvtkm