//============================================================================ // 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 2015 Sandia Corporation. // Copyright 2015 UT-Battelle, LLC. // Copyright 2015 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. //============================================================================ #ifndef vtk_m_cont_CellSetPermutation_h #define vtk_m_cont_CellSetPermutation_h #include #include #include #include #include #include #include #include #include #include #ifndef VTKM_DEFAULT_CELLSET_PERMUTATION_STORAGE_TAG #define VTKM_DEFAULT_CELLSET_PERMUTATION_STORAGE_TAG VTKM_DEFAULT_STORAGE_TAG #endif namespace vtkm { namespace cont { namespace internal { struct WriteConnectivity : public vtkm::worklet::WorkletMapPointToCell { typedef void ControlSignature(CellSetIn cellset, FieldInCell offset, WholeArrayOut<> connectivity); typedef void ExecutionSignature(PointCount, PointIndices, _2, _3); using InputDomain = _1; template 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 class PointToCell { private: using NumIndicesArrayType = decltype(make_ArrayHandlePermutation( std::declval().GetValidCellIds(), std::declval().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 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::ScanExclusive( vtkm::cont::make_ArrayHandleCast(conn.NumIndices, vtkm::Id()), conn.IndexOffsets); conn.Connectivity.Allocate(connectivityLength); vtkm::worklet::DispatcherMapTopology().Invoke( cellset, conn.IndexOffsets, conn.Connectivity); return conn; } }; // specialization for CellSetStructured template class PointToCell, PermutationArrayHandleType>> { private: using CellSetPermutationType = CellSetPermutation, PermutationArrayHandleType>; public: using ConnectivityArrays = vtkm::cont::internal::ConnectivityExplicitInternals< VTKM_DEFAULT_STORAGE_TAG, // dummy, shapes array is not used typename vtkm::cont::ArrayHandleConstant::StorageTag, VTKM_DEFAULT_STORAGE_TAG, typename vtkm::cont::ArrayHandleCounting::StorageTag>; template static ConnectivityArrays Get(const CellSetPermutationType& cellset, Device) { vtkm::Id numberOfCells = cellset.GetNumberOfCells(); vtkm::IdComponent numPointsInCell = vtkm::internal::ConnectivityStructuredInternals::NUM_POINTS_IN_CELL; vtkm::Id connectivityLength = numberOfCells * numPointsInCell; ConnectivityArrays conn; conn.NumIndices = make_ArrayHandleConstant(numPointsInCell, numberOfCells); conn.IndexOffsets = ArrayHandleCounting(0, numPointsInCell, numberOfCells); conn.Connectivity.Allocate(connectivityLength); vtkm::worklet::DispatcherMapTopology().Invoke( cellset, conn.IndexOffsets, conn.Connectivity); return conn; } }; template struct CellSetPermutationTraits; template struct CellSetPermutationTraits> { using OriginalCellSet = OriginalCellSet_; using PermutationArrayHandleType = PermutationArrayHandleType_; }; template struct CellSetPermutationTraits< CellSetPermutation, PermutationArrayHandleType_>> { using PreviousCellSet = CellSetPermutation; using PermutationArrayHandleType = vtkm::cont::ArrayHandlePermutation< PermutationArrayHandleType_, typename CellSetPermutationTraits::PermutationArrayHandleType>; using OriginalCellSet = typename CellSetPermutationTraits::OriginalCellSet; using Superclass = CellSetPermutation; }; } // internal template > class CellSetPermutation : public CellSet { public: using OriginalCellSetType = OriginalCellSetType_; using PermutationArrayHandleType = PermutationArrayHandleType_; VTKM_CONT CellSetPermutation(const PermutationArrayHandleType& validCellIds, const OriginalCellSetType& cellset, const std::string& name = std::string()) : CellSet(name) , ValidCellIds(validCellIds) , FullCellSet(cellset) { } VTKM_CONT CellSetPermutation(const std::string& name = std::string()) : CellSet(name) , ValidCellIds() , FullCellSet() { } VTKM_CONT const OriginalCellSetType& GetFullCellSet() const { return this->FullCellSet; } VTKM_CONT const PermutationArrayHandleType& GetValidCellIds() const { return this->ValidCellIds; } VTKM_CONT vtkm::Id GetNumberOfCells() const override { return this->ValidCellIds.GetNumberOfValues(); } 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 OriginalCellSetType& cellset) { this->ValidCellIds = validCellIds; this->FullCellSet = cellset; } VTKM_CONT vtkm::Id GetSchedulingRange(vtkm::TopologyElementTagCell) const { return this->ValidCellIds.GetNumberOfValues(); } VTKM_CONT vtkm::Id GetSchedulingRange(vtkm::TopologyElementTagPoint) const { return this->FullCellSet.GetNumberOfPoints(); } template struct ExecutionTypes; template struct ExecutionTypes { VTKM_IS_DEVICE_ADAPTER_TAG(Device); using ExecPortalType = typename PermutationArrayHandleType::template ExecutionTypes::PortalConst; using OrigExecObjectType = typename OriginalCellSetType::template ExecutionTypes< Device, vtkm::TopologyElementTagPoint, vtkm::TopologyElementTagCell>::ExecObjectType; using ExecObjectType = vtkm::exec::ConnectivityPermutedPointToCell; }; template struct ExecutionTypes { VTKM_IS_DEVICE_ADAPTER_TAG(Device); using ConnectiviyPortalType = typename vtkm::cont::ArrayHandle::template ExecutionTypes::PortalConst; using NumIndicesPortalType = typename vtkm::cont::ArrayHandle< vtkm::IdComponent>::template ExecutionTypes::PortalConst; using IndexOffsetPortalType = typename vtkm::cont::ArrayHandle::template ExecutionTypes::PortalConst; using ExecObjectType = vtkm::exec::ConnectivityPermutedCellToPoint; }; template VTKM_CONT typename ExecutionTypes::ExecObjectType PrepareForInput(Device device, vtkm::TopologyElementTagPoint from, vtkm::TopologyElementTagCell to) const { using ConnectivityType = typename ExecutionTypes::ExecObjectType; return ConnectivityType(this->ValidCellIds.PrepareForInput(device), this->FullCellSet.PrepareForInput(device, from, to)); } template VTKM_CONT typename ExecutionTypes::ExecObjectType PrepareForInput(Device device, vtkm::TopologyElementTagCell, vtkm::TopologyElementTagPoint) const { if (!this->CellToPoint.ElementsValid) { auto pointToCell = internal::PointToCell::Get(*this, device); internal::ComputeCellToPointConnectivity( this->CellToPoint, pointToCell, this->GetNumberOfPoints(), device); this->CellToPoint.BuildIndexOffsets(device); } using ConnectivityType = typename ExecutionTypes::ExecObjectType; return ConnectivityType(this->CellToPoint.Connectivity.PrepareForInput(device), this->CellToPoint.NumIndices.PrepareForInput(device), this->CellToPoint.IndexOffsets.PrepareForInput(device)); } VTKM_CONT void PrintSummary(std::ostream& out) const override { 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::StorageTag>; PermutationArrayHandleType ValidCellIds; OriginalCellSetType FullCellSet; mutable CellToPointConnectivity CellToPoint; }; template class CellSetPermutation, PermutationArrayHandleType2> : public internal::CellSetPermutationTraits< CellSetPermutation, PermutationArrayHandleType2>>::Superclass { private: using Superclass = typename internal::CellSetPermutationTraits::Superclass; public: VTKM_CONT CellSetPermutation(const PermutationArrayHandleType2& validCellIds, const CellSetPermutation& cellset, const std::string& name = std::string()) : Superclass(vtkm::cont::make_ArrayHandlePermutation(validCellIds, cellset.GetValidCellIds()), cellset.GetFullCellSet(), name) { } VTKM_CONT CellSetPermutation(const std::string& name = std::string()) : Superclass(name) { } VTKM_CONT void Fill(const PermutationArrayHandleType2& validCellIds, const CellSetPermutation& cellset) { this->ValidCellIds = make_ArrayHandlePermutation(validCellIds, cellset.GetValidCellIds()); this->FullCellSet = cellset.GetFullCellSet(); } }; template vtkm::cont::CellSetPermutation make_CellSetPermutation( const PermutationArrayHandleType& cellIndexMap, const OriginalCellSet& cellSet, const std::string& name) { VTKM_IS_CELL_SET(OriginalCellSet); VTKM_IS_ARRAY_HANDLE(PermutationArrayHandleType); return vtkm::cont::CellSetPermutation( cellIndexMap, cellSet, name); } template vtkm::cont::CellSetPermutation make_CellSetPermutation( const PermutationArrayHandleType& cellIndexMap, const OriginalCellSet& cellSet) { VTKM_IS_CELL_SET(OriginalCellSet); VTKM_IS_ARRAY_HANDLE(PermutationArrayHandleType); return vtkm::cont::make_CellSetPermutation(cellIndexMap, cellSet, cellSet.GetName()); } } } // namespace vtkm::cont #endif //vtk_m_cont_CellSetPermutation_h