Merge topic 'pointlocator'

ed3a64a5 Coding style improvment
7fa800b7 Update TestingPointLocatorUniformGrid.h
f1974cab Update TestingPointLocatorUniformGrid.h
508882fa PointLocatorUniformGrid

Acked-by: Kitware Robot <kwrobot@kitware.com>
Merge-request: !973
This commit is contained in:
Li-Ta Lo 2017-10-25 14:41:59 +00:00 committed by Kitware Robot
commit 3acd7c37a1
10 changed files with 468 additions and 0 deletions

@ -71,6 +71,7 @@ set(headers
Field.h
ImplicitFunction.h
MultiBlock.h
PointLocatorUniformGrid.h
RuntimeDeviceInformation.h
RuntimeDeviceTracker.h
Storage.h

@ -0,0 +1,237 @@
//============================================================================
// 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 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
// Copyright 2014 UT-Battelle, LLC.
// Copyright 2017 Los Alamos National Security.
//
// Under the terms of Contract DE-NA0003525 with NTESS,
// 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_PointLocatorUniformGrid_h
#define vtk_m_cont_PointLocatorUniformGrid_h
#include <vtkm/cont/ArrayHandleCounting.h>
#include <vtkm/cont/DeviceAdapter.h>
#include <vtkm/cont/DeviceAdapterAlgorithm.h>
#include <vtkm/worklet/DispatcherMapField.h>
#include <vtkm/worklet/WorkletMapField.h>
namespace vtkm
{
namespace worklet
{
template <typename T>
class PointLocatorUniformGrid
{
public:
PointLocatorUniformGrid(const vtkm::Vec<T, 3>& _min,
const vtkm::Vec<T, 3>& _max,
const vtkm::Vec<vtkm::Id, 3>& _dims)
: Min(_min)
, Max(_max)
, Dims(_dims)
{
}
class BinPointsWorklet : public vtkm::worklet::WorkletMapField
{
public:
typedef void ControlSignature(FieldIn<> coord, FieldOut<> label);
typedef void ExecutionSignature(_1, _2);
VTKM_CONT
BinPointsWorklet(vtkm::Vec<T, 3> _min, vtkm::Vec<T, 3> _max, vtkm::Vec<vtkm::Id, 3> _dims)
: Min(_min)
, Dims(_dims)
, Dxdydz((_max - Min) / Dims)
{
}
template <typename CoordVecType, typename IdType>
VTKM_EXEC void operator()(const CoordVecType& coord, IdType& label) const
{
vtkm::Vec<vtkm::Id, 3> ijk = (coord - Min) / Dxdydz;
label = ijk[0] + ijk[1] * Dims[0] + ijk[2] * Dims[0] * Dims[1];
}
private:
vtkm::Vec<T, 3> Min;
vtkm::Vec<vtkm::Id, 3> Dims;
vtkm::Vec<T, 3> Dxdydz;
};
class UniformGridSearch : public vtkm::worklet::WorkletMapField
{
public:
typedef void ControlSignature(FieldIn<> query,
WholeArrayIn<> coordIn,
WholeArrayIn<IdType> pointId,
WholeArrayIn<IdType> cellLower,
WholeArrayIn<IdType> cellUpper,
FieldOut<IdType> neighborId,
FieldOut<> distance);
typedef void ExecutionSignature(_1, _2, _3, _4, _5, _6, _7);
VTKM_CONT
UniformGridSearch(const vtkm::Vec<T, 3>& _min,
const vtkm::Vec<T, 3>& _max,
const vtkm::Vec<vtkm::Id, 3>& _dims)
: Min(_min)
, Dims(_dims)
, Dxdydz((_max - _min) / _dims)
{
}
template <typename CoordiVecType,
typename IdPortalType,
typename CoordiPortalType,
typename IdType,
typename CoordiType>
VTKM_EXEC void operator()(const CoordiVecType& queryCoord,
const CoordiPortalType& coordi_Handle,
const IdPortalType& pointId,
const IdPortalType& cellLower,
const IdPortalType& cellUpper,
IdType& nnId,
CoordiType& nnDis) const
{
auto nlayers = vtkm::Max(vtkm::Max(Dims[0], Dims[1]), Dims[2]);
vtkm::Vec<vtkm::Id, 3> xyz = (queryCoord - Min) / Dxdydz;
float min_distance = std::numeric_limits<float>::max();
vtkm::Id neareast = -1;
for (vtkm::Id layer = 0; layer < nlayers; layer++)
{
vtkm::Id minx = vtkm::Max(vtkm::Id(), xyz[0] - layer);
vtkm::Id maxx = vtkm::Min(Dims[0] - 1, xyz[0] + layer);
vtkm::Id miny = vtkm::Max(vtkm::Id(), xyz[1] - layer);
vtkm::Id maxy = vtkm::Min(Dims[1] - 1, xyz[1] + layer);
vtkm::Id minz = vtkm::Max(vtkm::Id(), xyz[2] - layer);
vtkm::Id maxz = vtkm::Min(Dims[2] - 1, xyz[2] + layer);
for (auto i = minx; i <= maxx; i++)
{
for (auto j = miny; j <= maxy; j++)
{
for (auto k = minz; k <= maxz; k++)
{
if (i == (xyz[0] + layer) || i == (xyz[0] - layer) || j == (xyz[1] + layer) ||
j == (xyz[1] - layer) || k == (xyz[2] + layer) || k == (xyz[2] - layer))
{
auto cellid = i + j * Dims[0] + k * Dims[0] * Dims[1];
auto lower = cellLower.Get(cellid);
auto upper = cellUpper.Get(cellid);
for (auto index = lower; index < upper; index++)
{
auto pointid = pointId.Get(index);
auto point = coordi_Handle.Get(pointid);
auto dx = point[0] - queryCoord[0];
auto dy = point[1] - queryCoord[1];
auto dz = point[2] - queryCoord[2];
auto distance2 = dx * dx + dy * dy + dz * dz;
if (distance2 < min_distance)
{
neareast = pointid;
min_distance = distance2;
nlayers = layer + 2;
}
}
}
}
}
}
}
nnId = neareast;
nnDis = vtkm::Sqrt(min_distance);
};
private:
vtkm::Vec<T, 3> Min;
vtkm::Vec<vtkm::Id, 3> Dims;
vtkm::Vec<T, 3> Dxdydz;
};
/// \brief Construct a 3D uniform grid for nearest neighbor search.
///
/// \param coords An ArrayHandle of x, y, z coordinates of input points.
/// \param device Tag for selecting device adapter
template <typename DeviceAdapter>
void Build(const vtkm::cont::ArrayHandle<vtkm::Vec<T, 3>>& coords, DeviceAdapter)
{
typedef vtkm::cont::DeviceAdapterAlgorithm<DeviceAdapter> Algorithm;
// generate unique id for each input point
vtkm::cont::ArrayHandleCounting<vtkm::Id> pointCounting(0, 1, coords.GetNumberOfValues());
Algorithm::Copy(pointCounting, PointIds);
// bin points into cells and give each of them the cell id.
BinPointsWorklet cellIdWorklet(Min, Max, Dims);
vtkm::worklet::DispatcherMapField<BinPointsWorklet> dispatchCellId(cellIdWorklet);
dispatchCellId.Invoke(coords, CellIds);
// Group points of the same cell together by sorting them according to the cell ids
Algorithm::SortByKey(CellIds, PointIds);
// for each cell, find the lower and upper bound of indices to the sorted point ids.
vtkm::cont::ArrayHandleCounting<vtkm::Id> cell_ids_counting(0, 1, Dims[0] * Dims[1] * Dims[2]);
Algorithm::UpperBounds(CellIds, cell_ids_counting, CellUpper);
Algorithm::LowerBounds(CellIds, cell_ids_counting, CellLower);
};
/// \brief Nearest neighbor search using a Uniform Grid
///
/// Parallel search of nearesat neighbor for each point in the \c queryPoints in the set of
/// \c coords. Returns neareast neighbot in \c nearestNeighborIds and distances to nearest
/// neighbor in \c distances.
///
/// \param coords Point coordinates for training dataset.
/// \param queryPoints Point coordinates to query for nearest neighbors.
/// \param nearestNeighborIds Neareast neighbor in the training dataset for each points in
/// the test set
/// \param distances Distance between query points and their nearest neighbors.
/// \param device Tag for selecting device adapter.
template <typename DeviceAdapter>
void FindNearestPoint(const vtkm::cont::ArrayHandle<vtkm::Vec<T, 3>>& coords,
const vtkm::cont::ArrayHandle<vtkm::Vec<T, 3>>& queryPoints,
vtkm::cont::ArrayHandle<vtkm::Id>& nearestNeighborIds,
vtkm::cont::ArrayHandle<T>& distances,
DeviceAdapter)
{
UniformGridSearch uniformGridSearch(Min, Max, Dims);
vtkm::worklet::DispatcherMapField<UniformGridSearch, DeviceAdapter> searchDispatcher(
uniformGridSearch);
searchDispatcher.Invoke(
queryPoints, coords, PointIds, CellLower, CellUpper, nearestNeighborIds, distances);
};
private:
vtkm::Vec<T, 3> Min;
vtkm::Vec<T, 3> Max;
vtkm::Vec<vtkm::Id, 3> Dims;
vtkm::cont::ArrayHandle<vtkm::Id> PointIds;
vtkm::cont::ArrayHandle<vtkm::Id> CellIds;
vtkm::cont::ArrayHandle<vtkm::Id> CellLower;
vtkm::cont::ArrayHandle<vtkm::Id> CellUpper;
};
}
}
#endif //vtk_m_cont_PointLocatorUniformGrid_h

@ -28,6 +28,7 @@ set(unit_tests
UnitTestCudaImplicitFunction.cu
UnitTestCudaMath.cu
UnitTestCudaShareUserProvidedManagedMemory.cu
UnitTestCudaPointLocatorUniformGrid.cxx
UnitTestCudaVirtualObjectCache.cu
)
vtkm_unit_tests(CUDA SOURCES ${unit_tests})

@ -0,0 +1,26 @@
//============================================================================
// 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 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
// Copyright 2017 UT-Battelle, LLC.
// Copyright 2017 Los Alamos National Security.
//
// Under the terms of Contract DE-NA0003525 with NTESS,
// 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/testing/TestingPointLocatorUniformGrid.h>
int UnitTestCudaPointLocatorUniformGrid(int, char* [])
{
return vtkm::cont::testing::Testing::Run(
TestingPointLocatorUniformGrid<VTKM_DEFAULT_DEVICE_ADAPTER_TAG>());
}

@ -27,6 +27,7 @@ set(unit_tests
UnitTestSerialDataSetSingleType.cxx
UnitTestSerialDeviceAdapter.cxx
UnitTestSerialImplicitFunction.cxx
UnitTestSerialPointLocatorUniformGrid.cxx
UnitTestSerialVirtualObjectCache.cxx
)
vtkm_unit_tests(TBB SOURCES ${unit_tests})

@ -0,0 +1,26 @@
//============================================================================
// 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 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
// Copyright 2017 UT-Battelle, LLC.
// Copyright 2017 Los Alamos National Security.
//
// Under the terms of Contract DE-NA0003525 with NTESS,
// 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/testing/TestingPointLocatorUniformGrid.h>
int UnitTestSerialPointLocatorUniformGrid(int, char* [])
{
return vtkm::cont::testing::Testing::Run(
TestingPointLocatorUniformGrid<VTKM_DEFAULT_DEVICE_ADAPTER_TAG>());
}

@ -27,6 +27,7 @@ set(unit_tests
UnitTestTBBDataSetSingleType.cxx
UnitTestTBBDeviceAdapter.cxx
UnitTestTBBImplicitFunction.cxx
UnitTestTBBPointLocatorUniformGrid.cxx
UnitTestTBBVirtualObjectCache.cxx
)
vtkm_unit_tests(TBB SOURCES ${unit_tests})

@ -0,0 +1,26 @@
//============================================================================
// 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 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
// Copyright 2017 UT-Battelle, LLC.
// Copyright 2017 Los Alamos National Security.
//
// Under the terms of Contract DE-NA0003525 with NTESS,
// 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/testing/TestingPointLocatorUniformGrid.h>
int UnitTestTBBPointLocatorUniformGrid(int, char* [])
{
return vtkm::cont::testing::Testing::Run(
TestingPointLocatorUniformGrid<VTKM_DEFAULT_DEVICE_ADAPTER_TAG>());
}

@ -30,6 +30,7 @@ set(headers
TestingDataSetSingleType.h
TestingFancyArrayHandles.h
TestingImplicitFunction.h
TestingPointLocatorUniformGrid.h
TestingVirtualObjectCache.h
)

@ -0,0 +1,148 @@
//============================================================================
// 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 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
// Copyright 2014 UT-Battelle, LLC.
// Copyright 2017 Los Alamos National Security.
//
// Under the terms of Contract DE-NA0003525 with NTESS,
// 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_testing_TestingPointLocatorUniformGrid_h
#define vtk_m_cont_testing_TestingPointLocatorUniformGrid_h
#include <random>
#include <vtkm/cont/testing/Testing.h>
#include <vtkm/cont/PointLocatorUniformGrid.h>
////brute force method /////
template <typename CoordiVecT, typename CoordiPortalT, typename CoordiT>
VTKM_EXEC_CONT vtkm::Id NNSVerify3D(CoordiVecT qc, CoordiPortalT coordiPortal, CoordiT& dis)
{
dis = std::numeric_limits<CoordiT>::max();
vtkm::Id nnpIdx = -1;
for (vtkm::Int32 i = 0; i < coordiPortal.GetNumberOfValues(); i++)
{
CoordiT splitX = coordiPortal.Get(i)[0];
CoordiT splitY = coordiPortal.Get(i)[1];
CoordiT splitZ = coordiPortal.Get(i)[2];
CoordiT _dis =
vtkm::Sqrt((splitX - qc[0]) * (splitX - qc[0]) + (splitY - qc[1]) * (splitY - qc[1]) +
(splitZ - qc[2]) * (splitZ - qc[2]));
if (_dis < dis)
{
dis = _dis;
nnpIdx = i;
}
}
return nnpIdx;
}
class NearestNeighborSearchBruteForce3DWorklet : public vtkm::worklet::WorkletMapField
{
public:
typedef void ControlSignature(FieldIn<> qcIn,
WholeArrayIn<> treeCoordiIn,
FieldOut<> nnIdOut,
FieldOut<> nnDisOut);
typedef void ExecutionSignature(_1, _2, _3, _4);
VTKM_CONT
NearestNeighborSearchBruteForce3DWorklet() {}
template <typename CoordiVecType, typename CoordiPortalType, typename IdType, typename CoordiType>
VTKM_EXEC void operator()(const CoordiVecType& qc,
const CoordiPortalType& coordiPortal,
IdType& nnId,
CoordiType& nnDis) const
{
nnDis = std::numeric_limits<CoordiType>::max();
nnId = NNSVerify3D(qc, coordiPortal, nnDis);
}
};
template <typename DeviceAdapter>
class TestingPointLocatorUniformGrid
{
public:
typedef vtkm::cont::DeviceAdapterAlgorithm<VTKM_DEFAULT_DEVICE_ADAPTER_TAG> Algorithm;
void TestTest() const
{
vtkm::Int32 nTrainingPoints = 1000;
vtkm::Int32 nTestingPoint = 1000;
std::vector<vtkm::Vec<vtkm::Float32, 3>> coordi;
///// randomly generate training points/////
std::default_random_engine dre;
std::uniform_real_distribution<vtkm::Float32> dr(0.0f, 10.0f);
for (vtkm::Int32 i = 0; i < nTrainingPoints; i++)
{
coordi.push_back(vtkm::make_Vec(dr(dre), dr(dre), dr(dre)));
}
auto coordi_Handle = vtkm::cont::make_ArrayHandle(coordi);
vtkm::worklet::PointLocatorUniformGrid<vtkm::Float32> uniformGrid(
{ 0.0f, 0.0f, 0.0f }, { 10.0f, 10.0f, 10.0f }, { 5, 5, 5 });
uniformGrid.Build(coordi_Handle, VTKM_DEFAULT_DEVICE_ADAPTER_TAG());
std::vector<vtkm::Vec<vtkm::Float32, 3>> qcVec;
for (vtkm::Int32 i = 0; i < nTestingPoint; i++)
{
qcVec.push_back(vtkm::make_Vec(dr(dre), dr(dre), dr(dre)));
}
auto qc_Handle = vtkm::cont::make_ArrayHandle(qcVec);
vtkm::cont::ArrayHandle<vtkm::Id> nnId_Handle;
vtkm::cont::ArrayHandle<vtkm::Float32> nnDis_Handle;
uniformGrid.FindNearestPoint(
coordi_Handle, qc_Handle, nnId_Handle, nnDis_Handle, VTKM_DEFAULT_DEVICE_ADAPTER_TAG());
vtkm::cont::ArrayHandle<vtkm::Id> bfnnId_Handle;
vtkm::cont::ArrayHandle<vtkm::Float32> bfnnDis_Handle;
NearestNeighborSearchBruteForce3DWorklet nnsbf3dWorklet;
vtkm::worklet::DispatcherMapField<NearestNeighborSearchBruteForce3DWorklet> nnsbf3DDispatcher(
nnsbf3dWorklet);
nnsbf3DDispatcher.Invoke(
qc_Handle, vtkm::cont::make_ArrayHandle(coordi), bfnnId_Handle, bfnnDis_Handle);
///// verfity search result /////
bool passTest = true;
for (vtkm::Int32 i = 0; i < nTestingPoint; i++)
{
vtkm::Id workletIdx = nnId_Handle.GetPortalControl().Get(i);
vtkm::Float32 workletDis = nnDis_Handle.GetPortalConstControl().Get(i);
vtkm::Id bfworkletIdx = bfnnId_Handle.GetPortalControl().Get(i);
vtkm::Float32 bfworkletDis = bfnnDis_Handle.GetPortalConstControl().Get(i);
if (workletIdx != bfworkletIdx)
{
std::cout << "bf index: " << bfworkletIdx << ", dis: " << bfworkletDis
<< ", grid: " << workletIdx << ", dis " << workletDis << std::endl;
passTest = false;
}
}
VTKM_TEST_ASSERT(passTest, "Uniform Grid NN search result incorrect.");
}
void operator()() const { this->TestTest(); }
};
#endif // vtk_m_cont_testing_TestingPointLocatorUniformGrid_h