//
// ---------------------------------------------------------------------
-#ifndef _dealii__numerics_kdtree_distance_h
-#define _dealii__numerics_kdtree_distance_h
+#ifndef _dealii__numerics_kdtree_h
+#define _dealii__numerics_kdtree_h
#include <deal.II/base/config.h>
* A wrapper for the nanoflann library, used to compute the distance from a
* collection of points, and to efficiently return nearest neighbors to a
* target point. This function uses nanoflann to efficiently partition the
- * space in a tree. The cost of each query is then roughly of order log(n),
- * where n is the number of points stored in this class.
+ * space in a k-dimensional tree. The cost of each query is then roughly of
+ * order log(n), where n is the number of points stored in this class.
*
- * The wrapper provides methods that give access to some of thefunctionalities
+ * The wrapper provides methods that give access to some of the functionalities
* of the nanoflann library, like searching the n nearest neighbors, or
- * searching the points that fall within a raidius of a target point.
+ * searching the points that fall within a radius of a target point.
+ *
+ * @quotation
+ * From wikipedia (https://en.wikipedia.org/wiki/K-d_tree):
+ *
+ * A k-d tree is a binary tree in which every node is a k-dimensional point.
+ * Every non-leaf node can be thought of as implicitly generating a splitting
+ * hyperplane that divides the space into two parts, known as half-spaces.
+ * Points to the left of this hyperplane are represented by the left subtree of
+ * that node and points right of the hyperplane are represented by the right
+ * subtree. The hyperplane direction is chosen in the following way: every node
+ * in the tree is associated with one of the k-dimensions, with the hyperplane
+ * perpendicular to that dimension's axis. So, for example, if for a particular
+ * split the "x" axis is chosen, all points in the subtree with a smaller "x"
+ * value than the node will appear in the left subtree and all points with
+ * larger "x" value will be in the right subtree. In such a case, the
+ * hyperplane would be set by the x-value of the point, and its normal would be
+ * the unit x-axis.
+ * @endquotation
+ *
+ * @author Luca Heltai, 2017.
*/
template<int dim>
-class KDTreeDistance
+class KDTree
{
public:
/**
* a vector of points will result in an exception. Only a reference
* to the points is stored, so you should make sure that the life of
* the the vector you pass is longer than the life of this class, or
- * you'll get undefinite behaviour.
+ * you'll get undefined behaviour.
*
* If you update the vector of points in someway, remember to call
* again the set_points() method. The tree and the index are
* your points, and do not call again set_points(), your results
* will likely be wrong.
*/
- KDTreeDistance(const unsigned int &max_leaf_size=10,
- const std::vector<Point<dim> > &pts=std::vector<Point<dim> >());
+ KDTree(const unsigned int &max_leaf_size=10,
+ const std::vector<Point<dim> > &pts=std::vector<Point<dim> >());
/**
* Reference to the vector of points from which we want to compute
* the distance.
*/
- const std::vector<Point<dim> > &points; //!< A const ref to the data set origin
+ const std::vector<Point<dim> > &points;
/**
- * The constrcutor needs the data set source.
+ * The constructor needs the vector of points from which we want to build
+ * the tree.
*/
- PointCloudAdaptor(const std::vector<Point<dim> > &_points) : points(_points) { }
+ PointCloudAdaptor(const std::vector<Point<dim> > &_points);
/**
* Return number of points in the data set (required by nanoflann).
*/
- inline size_t kdtree_get_point_count() const
- {
- return points.size();
- }
+ size_t kdtree_get_point_count() const;
/**
* Return the L2 distance between points
*/
- inline coord_t kdtree_distance(const coord_t *p1, const size_t idx_p2,size_t size) const
- {
- AssertDimension(size, dim);
- coord_t res=0.0;
- for (size_t d=0; d<size; ++d)
- res += (p1[d]-points[idx_p2][d])*(p1[d]-points[idx_p2][d]);
- return std::sqrt(res);
- }
+ coord_t kdtree_distance(const coord_t *p1, const size_t idx_p2,size_t size) const;
/**
- * Return the dim'th component of the idx'th point in the class.
+ * Return the d-th component of the idx-th point in the class.
*/
- inline coord_t kdtree_get_pt(const size_t idx, int d) const
- {
- AssertIndexRange(d,dim);
- return points[idx][d];
- }
+ coord_t kdtree_get_pt(const size_t idx, int d) const;
/**
* A typedef for the actual KDTree object.
*/
typedef typename nanoflann::KDTreeSingleIndexAdaptor<nanoflann::L2_Simple_Adaptor<double, PointCloudAdaptor> ,
- PointCloudAdaptor, dim, unsigned int> KDTree;
+ PointCloudAdaptor, dim, unsigned int> NanoFlannKDTree;
/**
* @param[in] point: the target point
* @param[in] radius: the radius of the ball
* @param[in] sorted: sort the output results in ascending order with respect to distances
- * @param[out] matches:
*
* @return vector of indices and distances of the matching points
*/
bool sorted=false) const;
/**
- * Fill two vectors with the indices and distances of the closest
- * points to the given target point. The vectors are filled with
- * indices and distances until there is space in them. You should
- * resize them to the number of closest points you wish to get. An
- * assertion is thrown if the vectors do not have the same size.
+ * Fill a vectors with the indices and distances of the closest @p n_points
+ * points to the given target point.
*
* @param[in] target: the target point
* @param[in] n_points: the number of requested points
/**
* The actual kdtree.
*/
- std::unique_ptr<KDTree> kdtree;
+ std::unique_ptr<NanoFlannKDTree> kdtree;
};
//------------ inline functions -------------
+#ifndef DOXYGEN
template<int dim>
inline
-unsigned int KDTreeDistance<dim>::size() const
+unsigned int KDTree<dim>::size() const
{
if (adaptor)
return adaptor->points.size();
return 0;
};
+
+
template<int dim>
inline const Point<dim> &
-KDTreeDistance<dim>::operator[](unsigned int i) const
+KDTree<dim>::operator[](unsigned int i) const
{
AssertIndexRange(i, size());
return adaptor->points[i];
}
+
+template<int dim>
+KDTree<dim>::PointCloudAdaptor::PointCloudAdaptor(const std::vector<Point<dim> > &_points)
+ : points(_points)
+{}
+
+
+
+template<int dim>
+inline size_t
+KDTree<dim>::PointCloudAdaptor::kdtree_get_point_count() const
+{
+ return points.size();
+}
+
+
+
+template<int dim>
+inline double
+KDTree<dim>::PointCloudAdaptor::kdtree_get_pt(const size_t idx, int d) const
+{
+ AssertIndexRange(d,dim);
+ return points[idx][d];
+}
+
+
+
template<int dim>
template <class BBOX>
inline bool
-KDTreeDistance<dim>::PointCloudAdaptor::kdtree_get_bbox(BBOX &) const
+KDTree<dim>::PointCloudAdaptor::kdtree_get_bbox(BBOX &) const
{
return false;
}
+
+
+template<int dim>
+inline double
+KDTree<dim>::PointCloudAdaptor::kdtree_distance(const double *p1, const size_t idx_p2,size_t size) const
+{
+ AssertDimension(size, dim);
+ double res=0.0;
+ for (size_t d=0; d<size; ++d)
+ res += (p1[d]-points[idx_p2][d])*(p1[d]-points[idx_p2][d]);
+ return std::sqrt(res);
+}
+#endif
+
DEAL_II_NAMESPACE_CLOSE
# endif // DEAL_II_WITH_NANO_FLANN
-#include <deal.II/numerics/kdtree_distance.h>
+#include <deal.II/numerics/kdtree.h>
#ifdef DEAL_II_WITH_NANOFLANN
template<int dim>
-KDTreeDistance<dim>::KDTreeDistance(const unsigned int &max_leaf_size,
- const std::vector<Point<dim> > &pts)
+KDTree<dim>::KDTree(const unsigned int &max_leaf_size,
+ const std::vector<Point<dim> > &pts)
: max_leaf_size(max_leaf_size)
{
if (pts.size() > 0)
template<int dim>
-std::vector<std::pair<unsigned int, double> > KDTreeDistance<dim>::get_points_within_ball(const Point<dim> ¢er,
+std::vector<std::pair<unsigned int, double> > KDTree<dim>::get_points_within_ball(const Point<dim> ¢er,
const double &radius,
bool sorted) const
{
}
template<int dim>
-std::vector<std::pair<unsigned int, double> > KDTreeDistance<dim>::get_closest_points(const Point<dim> &target,
+std::vector<std::pair<unsigned int, double> > KDTree<dim>::get_closest_points(const Point<dim> &target,
const unsigned int n_points) const
{
Assert(adaptor, ExcNotInitialized());
}
template<int dim>
-void KDTreeDistance<dim>::set_points(const std::vector<Point<dim> > &pts)
+void KDTree<dim>::set_points(const std::vector<Point<dim> > &pts)
{
Assert(pts.size() > 0, ExcMessage("Expecting a non zero set of points."));
adaptor = std_cxx14::make_unique<PointCloudAdaptor>(pts);
- kdtree = std_cxx14::make_unique<KDTree>(dim, *adaptor, nanoflann::KDTreeSingleIndexAdaptorParams(max_leaf_size));
+ kdtree = std_cxx14::make_unique<NanoFlannKDTree>(dim, *adaptor, nanoflann::KDTreeSingleIndexAdaptorParams(max_leaf_size));
kdtree->buildIndex();
}
-template class KDTreeDistance<1>;
-template class KDTreeDistance<2>;
-template class KDTreeDistance<3>;
+template class KDTree<1>;
+template class KDTree<2>;
+template class KDTree<3>;
DEAL_II_NAMESPACE_CLOSE
-//-----------------------------------------------------------
+// ---------------------------------------------------------------------
//
-// Copyright (C) 2015 by the deal2lkit authors
+// Copyright (C) 2017 by the deal.II authors
//
-// This file is part of the deal2lkit library.
+// This file is part of the deal.II library.
//
-// The deal2lkit library is free software; you can use it, redistribute
-// it, and/or modify it under the terms of the GNU Lesser General
-// Public License as published by the Free Software Foundation; either
-// version 2.1 of the License, or (at your option) any later version.
-// The full text of the license can be found in the file LICENSE at
-// the top level of the deal2lkit distribution.
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal2lkit distribution.
//
-//-----------------------------------------------------------
+//---------------------------------------------------------------------
// Create a list of points, and compute the minimum distance from some other points
// to this set, using a kdtree library
#include "../tests.h"
-#include <deal.II/numerics/kdtree_distance.h>
+#include <deal.II/numerics/kdtree.h>
int main ()
-//-----------------------------------------------------------
+// ---------------------------------------------------------------------
//
-// Copyright (C) 2015 by the deal2lkit authors
+// Copyright (C) 2017 by the deal.II authors
//
-// This file is part of the deal2lkit library.
+// This file is part of the deal.II library.
//
-// The deal2lkit library is free software; you can use it, redistribute
-// it, and/or modify it under the terms of the GNU Lesser General
-// Public License as published by the Free Software Foundation; either
-// version 2.1 of the License, or (at your option) any later version.
-// The full text of the license can be found in the file LICENSE at
-// the top level of the deal2lkit distribution.
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal2lkit distribution.
//
-//-----------------------------------------------------------
+//---------------------------------------------------------------------
// Create a list of points, and compute the closest points to a given one.
#include "../tests.h"
-#include <deal.II/numerics/kdtree_distance.h>
+#include <deal.II/numerics/kdtree.h>
int main ()
{
-//-----------------------------------------------------------
+// ---------------------------------------------------------------------
//
-// Copyright (C) 2015 by the deal2lkit authors
+// Copyright (C) 2017 by the deal.II authors
//
-// This file is part of the deal2lkit library.
+// This file is part of the deal.II library.
//
-// The deal2lkit library is free software; you can use it, redistribute
-// it, and/or modify it under the terms of the GNU Lesser General
-// Public License as published by the Free Software Foundation; either
-// version 2.1 of the License, or (at your option) any later version.
-// The full text of the license can be found in the file LICENSE at
-// the top level of the deal2lkit distribution.
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal2lkit distribution.
//
-//-----------------------------------------------------------
+//---------------------------------------------------------------------
// Create a list of points, and the ones within a specified radius of the target points
#include "../tests.h"
-#include <deal.II/numerics/kdtree_distance.h>
+#include <deal.II/numerics/kdtree.h>
int main ()
{