]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Addressed comments, renamed classes, externalized implementation
authorLuca Heltai <luca.heltai@sissa.it>
Mon, 14 Aug 2017 00:32:10 +0000 (18:32 -0600)
committerLuca Heltai <luca.heltai@sissa.it>
Mon, 14 Aug 2017 00:32:10 +0000 (18:32 -0600)
12 files changed:
cmake/configure/configure_nanoflann.cmake
cmake/modules/FindNANOFLANN.cmake
doc/news/changes/major/20170813LucaHeltai
include/deal.II/numerics/kdtree.h [moved from include/deal.II/numerics/kdtree_distance.h with 65% similarity]
source/numerics/CMakeLists.txt
source/numerics/kdtree.cc [moved from source/numerics/kdtree_distance.cc with 65% similarity]
tests/numerics/kdtree_01.cc [moved from tests/numerics/kdtree_distance_01.cc with 62% similarity]
tests/numerics/kdtree_01.output [moved from tests/numerics/kdtree_distance_01.output with 100% similarity]
tests/numerics/kdtree_02.cc [moved from tests/numerics/kdtree_distance_02.cc with 63% similarity]
tests/numerics/kdtree_02.output [moved from tests/numerics/kdtree_distance_02.output with 100% similarity]
tests/numerics/kdtree_03.cc [moved from tests/numerics/kdtree_distance_03.cc with 63% similarity]
tests/numerics/kdtree_03.output [moved from tests/numerics/kdtree_distance_03.output with 100% similarity]

index c6bc02b6934aa8ab384bcb91ab42785e293c03fd..6585109417d7aa65e3722c7e9413b4165e844ffd 100644 (file)
@@ -1,6 +1,6 @@
 ## ---------------------------------------------------------------------
 ##
-## Copyright (C) 2012 - 2014 by the deal.II authors
+## Copyright (C) 2017 by the deal.II authors
 ##
 ## This file is part of the deal.II library.
 ##
index b34d80a11d295397df2143110fa8ff2f8f416030..79bb99c089d1a9d888c3b7d1cbec00b6802460a7 100644 (file)
@@ -1,6 +1,6 @@
 ## ---------------------------------------------------------------------
 ##
-## Copyright (C) 2014 - 2015 by the deal.II authors
+## Copyright (C) 2017 by the deal.II authors
 ##
 ## This file is part of the deal.II library.
 ##
index a1ccfab047fa624696b65503348eecd142534549..ccccf9c5340a1b38a6fe45a6175bac3b5e797902 100644 (file)
@@ -1,7 +1,7 @@
 New: A new KDTreeDistance class has been added that interfaces
-to the nanoflann library (https://github.com/jlblancoc/nanoflann).
+the nanoflann library (https://github.com/jlblancoc/nanoflann).
 This class can be used to extract nearest neighbour information 
-on collection of points, query for the closes points to a target
+on collection of points, query for the closest points to a target
 point or all points contained within a given distance.
 <br>
 (Luca Heltai, 2017/08/13)
similarity index 65%
rename from include/deal.II/numerics/kdtree_distance.h
rename to include/deal.II/numerics/kdtree.h
index 08e048aa103ae7a8365306c015172dbf3f3e6c76..9e9556375a77ef5d749bf21d5138d55e5dea09c0 100644 (file)
@@ -13,8 +13,8 @@
 //
 // ---------------------------------------------------------------------
 
-#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>
 
@@ -33,15 +33,35 @@ DEAL_II_NAMESPACE_OPEN
  * 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:
   /**
@@ -56,7 +76,7 @@ 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
@@ -65,8 +85,8 @@ public:
    * 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> >());
 
 
   /**
@@ -86,45 +106,32 @@ public:
      * 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;
 
 
     /**
@@ -143,7 +150,7 @@ public:
    * 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;
 
 
   /**
@@ -188,7 +195,6 @@ public:
    * @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
    */
@@ -197,11 +203,8 @@ public:
       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
@@ -227,15 +230,16 @@ private:
   /**
    * 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();
@@ -243,23 +247,66 @@ unsigned int KDTreeDistance<dim>::size() const
     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
index 0d5252104a318112bc550acbb89f9cb9ac4125fc..fd906c011a9908116d5330eb2805b98fd9fb70e4 100644 (file)
@@ -23,7 +23,7 @@ SET(_unity_include_src
   data_postprocessor.cc
   dof_output_operator.cc
   histogram.cc
-  kdtree_distance.cc
+  kdtree.cc
   matrix_tools_once.cc
   matrix_tools.cc
   time_dependent.cc
similarity index 65%
rename from source/numerics/kdtree_distance.cc
rename to source/numerics/kdtree.cc
index 04345d5125f2e3a5c94a13806a81a959a9f7a936..39676cd65433043838d8d17b82b340db0fcca1d8 100644 (file)
@@ -1,4 +1,4 @@
-#include <deal.II/numerics/kdtree_distance.h>
+#include <deal.II/numerics/kdtree.h>
 
 #ifdef DEAL_II_WITH_NANOFLANN
 
@@ -8,8 +8,8 @@ DEAL_II_NAMESPACE_OPEN
 
 
 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)
@@ -18,7 +18,7 @@ KDTreeDistance<dim>::KDTreeDistance(const unsigned int &max_leaf_size,
 
 
 template<int dim>
-std::vector<std::pair<unsigned int, double> > KDTreeDistance<dim>::get_points_within_ball(const Point<dim> &center,
+std::vector<std::pair<unsigned int, double> > KDTree<dim>::get_points_within_ball(const Point<dim> &center,
     const double &radius,
     bool sorted) const
 {
@@ -36,7 +36,7 @@ std::vector<std::pair<unsigned int, double> > KDTreeDistance<dim>::get_points_wi
 }
 
 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());
@@ -52,18 +52,18 @@ std::vector<std::pair<unsigned int, double> > KDTreeDistance<dim>::get_closest_p
 }
 
 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
 
similarity index 62%
rename from tests/numerics/kdtree_distance_01.cc
rename to tests/numerics/kdtree_01.cc
index 2c0ae731e67fa166de612a1dd4e9d5083f1dadc0..9879dd38a6ec94cc13ffcf8c3cf8dcfc33ab2a8a 100644 (file)
@@ -1,23 +1,23 @@
-//-----------------------------------------------------------
+// ---------------------------------------------------------------------
 //
-//    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 ()
similarity index 63%
rename from tests/numerics/kdtree_distance_02.cc
rename to tests/numerics/kdtree_02.cc
index 0706c20825a103638edab5fa75fb474a3e175801..ff423b7f73a058a84416622633cc24f5ccfe6118 100644 (file)
@@ -1,22 +1,22 @@
-//-----------------------------------------------------------
+// ---------------------------------------------------------------------
 //
-//    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 ()
 {
similarity index 63%
rename from tests/numerics/kdtree_distance_03.cc
rename to tests/numerics/kdtree_03.cc
index df87f83f2e8f939c6a2f548f246d33445f86a7b5..b7fefa610de11522e3318920bc5158c8d61633b4 100644 (file)
@@ -1,22 +1,22 @@
-//-----------------------------------------------------------
+// ---------------------------------------------------------------------
 //
-//    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 ()
 {

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.