/**
- * Adaptor class used internally by nanoflann. Class actually stores
- * a reference to the vector of points, and generates some helper
- * functions for nanoflann.
+ * Adaptor class used internally by nanoflann. This class stores a reference
+ * to the vector of points, and generates some helper functions for
+ * nanoflann.
*/
struct PointCloudAdaptor
{
*
* @param[in] point: the target point
* @param[in] radius: the radius of the ball
- * @param[out] mathes: indices and distances of the matching points
* @param[in] sorted: sort the output results in ascending order with respect to distances
+ * @param[out] matches:
*
- * @return number of points that are within the ball
+ * @return vector of indices and distances of the matching points
*/
- unsigned int get_points_within_ball(const Point<dim> &target, const double &radius,
- std::vector<std::pair<unsigned int, double> > &matches,
- bool sorted=false) const;
+ std::vector<std::pair<unsigned int, double> > get_points_within_ball(const Point<dim> &target,
+ const double &radius,
+ bool sorted=false) const;
/**
* Fill two vectors with the indices and distances of the closest
* assertion is thrown if the vectors do not have the same size.
*
* @param[in] target: the target point
- * @param[out] indices: indices of the matching points
- * @param[out] distances: distances of the matching points
+ * @param[in] n_points: the number of requested points
+ *
+ * @return a vector of pairs of indices and distances of the matching points
*/
- void get_closest_points(const Point<dim> &target,
- std::vector<unsigned int> &indices,
- std::vector<double> &distances) const;
+ std::vector<std::pair<unsigned int, double> > get_closest_points(const Point<dim> &target,
+ const unsigned int n_points) const;
private:
/**
template<int dim>
-unsigned int KDTreeDistance<dim>::get_points_within_ball(const Point<dim> ¢er, const double &radius,
- std::vector<std::pair<unsigned int, double> > &matches,
- bool sorted) const
+std::vector<std::pair<unsigned int, double> > KDTreeDistance<dim>::get_points_within_ball(const Point<dim> ¢er,
+ const double &radius,
+ bool sorted) const
{
+ std::vector<std::pair<unsigned int, double> > matches;
Assert(adaptor, ExcNotInitialized());
Assert(kdtree, ExcInternalError());
nanoflann::SearchParams params;
params.sorted = sorted;
- return kdtree->radiusSearch(¢er[0], radius, matches, params);
+ kdtree->radiusSearch(¢er[0], radius, matches, params);
+ return matches;
}
template<int dim>
-void KDTreeDistance<dim>::get_closest_points(const Point<dim> &target,
- std::vector<unsigned int> &indices,
- std::vector<double> &distances) const
+std::vector<std::pair<unsigned int, double> > KDTreeDistance<dim>::get_closest_points(const Point<dim> &target,
+ const unsigned int n_points) const
{
Assert(adaptor, ExcNotInitialized());
Assert(kdtree, ExcInternalError());
- AssertDimension(indices.size(), distances.size());
-
- kdtree->knnSearch(&target[0], indices.size(), &indices[0], &distances[0]);
+ std::vector<unsigned int> indices(n_points);
+ std::vector<double> distances(n_points);
+ std::vector<std::pair<unsigned int, double> > matches(n_points);
+
+ kdtree->knnSearch(&target[0], n_points, &indices[0], &distances[0]);
+ for(unsigned int i=0; i<n_points; ++i)
+ matches[i] = std::make_pair(indices[i], distances[i]);
+ return matches;
}
template<int dim>
--- /dev/null
+//-----------------------------------------------------------
+//
+// Copyright (C) 2015 by the deal2lkit authors
+//
+// This file is part of the deal2lkit 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.
+//
+//-----------------------------------------------------------
+
+// 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>
+
+
+int main ()
+{
+ initlog();
+
+ KDTreeDistance<2> kdtree;
+
+ std::vector<Point<2> > points;
+
+ // Add four points
+ points.push_back(Point<2>(0,0));
+ points.push_back(Point<2>(0,1));
+ points.push_back(Point<2>(1,0));
+ points.push_back(Point<2>(1,1));
+
+ deallog << "Distance from unit square:" << std::endl;
+
+ std::vector<Point<2> > test_points;
+ test_points.push_back(Point<2>(.5, .5));
+ test_points.push_back(Point<2>(2, 0));
+ test_points.push_back(Point<2>(2, 2));
+
+ kdtree.set_points(points);
+
+ for (auto &p : test_points) {
+ auto res = kdtree.get_closest_points(p,1)[0];
+ deallog << "P: " << p << ", distance: " << res.second << ", index: " << res.first << std::endl;
+ }
+
+ deallog << "Consistency checking: the following are all the points in the set." << std::endl;
+ for (auto &p : points) {
+ auto res = kdtree.get_closest_points(p,1)[0];
+ deallog << "P: " << p << ", distance: " << res.second << ", index: " << res.first << std::endl;
+ }
+}
--- /dev/null
+
+DEAL::Distance from unit square:
+DEAL::P: 0.500000 0.500000, distance: 0.707107, index: 0
+DEAL::P: 2.00000 0.00000, distance: 1.00000, index: 2
+DEAL::P: 2.00000 2.00000, distance: 1.41421, index: 3
+DEAL::Consistency checking: the following are all the points in the set.
+DEAL::P: 0.00000 0.00000, distance: 0.00000, index: 0
+DEAL::P: 0.00000 1.00000, distance: 0.00000, index: 1
+DEAL::P: 1.00000 0.00000, distance: 0.00000, index: 2
+DEAL::P: 1.00000 1.00000, distance: 0.00000, index: 3
--- /dev/null
+//-----------------------------------------------------------
+//
+// Copyright (C) 2015 by the deal2lkit authors
+//
+// This file is part of the deal2lkit 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.
+//
+//-----------------------------------------------------------
+
+// Create a list of points, and compute the closest points to a given one.
+
+#include "../tests.h"
+#include <deal.II/numerics/kdtree_distance.h>
+
+int main ()
+{
+ initlog();
+
+ KDTreeDistance<2> kdtree;
+
+ std::vector<Point<2> > points;
+
+ // Add four points
+ points.push_back(Point<2>(0,0));
+ points.push_back(Point<2>(0,1));
+ points.push_back(Point<2>(1,0));
+ points.push_back(Point<2>(1,1));
+
+ std::vector<Point<2> > test_points;
+ test_points.push_back(Point<2>(.5, .5));
+ test_points.push_back(Point<2>(2, 0));
+ test_points.push_back(Point<2>(2, 2));
+
+ kdtree.set_points(points);
+
+ std::vector<unsigned int> indices;
+ std::vector<double> distances;
+
+ // Get closest points. Do a few rounds
+ for (auto &p : test_points)
+ {
+ for (unsigned int i=1; i<points.size()+1; ++i)
+ {
+ auto res = kdtree.get_closest_points(p, i);
+
+ deallog << std::endl << "The first " << i << " closest points to " << p << " are:" << std::endl;
+ for (unsigned int j=0; j<i; ++j)
+ {
+ deallog << "points[" << res[j].first << "] = " << points[res[j].first] << ", distance: "
+ << res[j].second << std::endl;
+ }
+ }
+ }
+}
--- /dev/null
+
+DEAL::
+DEAL::The first 1 closest points to 0.500000 0.500000 are:
+DEAL::points[0] = 0.00000 0.00000, distance: 0.707107
+DEAL::
+DEAL::The first 2 closest points to 0.500000 0.500000 are:
+DEAL::points[0] = 0.00000 0.00000, distance: 0.707107
+DEAL::points[1] = 0.00000 1.00000, distance: 0.707107
+DEAL::
+DEAL::The first 3 closest points to 0.500000 0.500000 are:
+DEAL::points[0] = 0.00000 0.00000, distance: 0.707107
+DEAL::points[1] = 0.00000 1.00000, distance: 0.707107
+DEAL::points[2] = 1.00000 0.00000, distance: 0.707107
+DEAL::
+DEAL::The first 4 closest points to 0.500000 0.500000 are:
+DEAL::points[0] = 0.00000 0.00000, distance: 0.707107
+DEAL::points[1] = 0.00000 1.00000, distance: 0.707107
+DEAL::points[2] = 1.00000 0.00000, distance: 0.707107
+DEAL::points[3] = 1.00000 1.00000, distance: 0.707107
+DEAL::
+DEAL::The first 1 closest points to 2.00000 0.00000 are:
+DEAL::points[2] = 1.00000 0.00000, distance: 1.00000
+DEAL::
+DEAL::The first 2 closest points to 2.00000 0.00000 are:
+DEAL::points[2] = 1.00000 0.00000, distance: 1.00000
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::
+DEAL::The first 3 closest points to 2.00000 0.00000 are:
+DEAL::points[2] = 1.00000 0.00000, distance: 1.00000
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::points[0] = 0.00000 0.00000, distance: 2.00000
+DEAL::
+DEAL::The first 4 closest points to 2.00000 0.00000 are:
+DEAL::points[2] = 1.00000 0.00000, distance: 1.00000
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::points[0] = 0.00000 0.00000, distance: 2.00000
+DEAL::points[1] = 0.00000 1.00000, distance: 2.23607
+DEAL::
+DEAL::The first 1 closest points to 2.00000 2.00000 are:
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::
+DEAL::The first 2 closest points to 2.00000 2.00000 are:
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::points[1] = 0.00000 1.00000, distance: 2.23607
+DEAL::
+DEAL::The first 3 closest points to 2.00000 2.00000 are:
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::points[1] = 0.00000 1.00000, distance: 2.23607
+DEAL::points[2] = 1.00000 0.00000, distance: 2.23607
+DEAL::
+DEAL::The first 4 closest points to 2.00000 2.00000 are:
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::points[1] = 0.00000 1.00000, distance: 2.23607
+DEAL::points[2] = 1.00000 0.00000, distance: 2.23607
+DEAL::points[0] = 0.00000 0.00000, distance: 2.82843
--- /dev/null
+//-----------------------------------------------------------
+//
+// Copyright (C) 2015 by the deal2lkit authors
+//
+// This file is part of the deal2lkit 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.
+//
+//-----------------------------------------------------------
+
+// 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>
+
+int main ()
+{
+ initlog();
+
+ KDTreeDistance<2> kdtree;
+
+ std::vector<Point<2> > points;
+
+ // Add four points
+ points.push_back(Point<2>(0,0));
+ points.push_back(Point<2>(0,1));
+ points.push_back(Point<2>(1,0));
+ points.push_back(Point<2>(1,1));
+
+ std::vector<Point<2> > test_points;
+ test_points.push_back(Point<2>(2, 0));
+ test_points.push_back(Point<2>(2, 2));
+
+ std::vector<double> radii;
+ radii.push_back(.8);
+ radii.push_back(1.1);
+ radii.push_back(1.5);
+ radii.push_back(2);
+
+
+ kdtree.set_points(points);
+
+ std::vector<std::pair<unsigned int, double> > matches;
+
+ // Get points within ball
+ for (auto &p : test_points)
+ {
+ for (auto &r : radii)
+ {
+
+ auto matches = kdtree.get_points_within_ball(p, r, true);
+
+ deallog << std::endl << "At distance less than " << r << " from " << p << " we have:" << std::endl;
+ for (auto &m : matches)
+ {
+ deallog << "points[" << m.first << "] = " << points[m.first] << ", distance: "
+ << m.second << std::endl;
+ }
+ }
+ }
+}
--- /dev/null
+
+DEAL::
+DEAL::At distance less than 0.800000 from 2.00000 0.00000 we have:
+DEAL::
+DEAL::At distance less than 1.10000 from 2.00000 0.00000 we have:
+DEAL::points[2] = 1.00000 0.00000, distance: 1.00000
+DEAL::
+DEAL::At distance less than 1.50000 from 2.00000 0.00000 we have:
+DEAL::points[2] = 1.00000 0.00000, distance: 1.00000
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::
+DEAL::At distance less than 2.00000 from 2.00000 0.00000 we have:
+DEAL::points[2] = 1.00000 0.00000, distance: 1.00000
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::
+DEAL::At distance less than 0.800000 from 2.00000 2.00000 we have:
+DEAL::
+DEAL::At distance less than 1.10000 from 2.00000 2.00000 we have:
+DEAL::
+DEAL::At distance less than 1.50000 from 2.00000 2.00000 we have:
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421
+DEAL::
+DEAL::At distance less than 2.00000 from 2.00000 2.00000 we have:
+DEAL::points[3] = 1.00000 1.00000, distance: 1.41421