From: Luca Heltai Date: Sun, 13 Aug 2017 23:30:22 +0000 (-0600) Subject: Changed output arguments and added 3 tests. X-Git-Tag: v9.0.0-rc1~1262^2~3 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=ba3d4e65afd79d49c615092050751959009c4227;p=dealii.git Changed output arguments and added 3 tests. --- diff --git a/include/deal.II/numerics/kdtree_distance.h b/include/deal.II/numerics/kdtree_distance.h index 810018eafa..9676796531 100644 --- a/include/deal.II/numerics/kdtree_distance.h +++ b/include/deal.II/numerics/kdtree_distance.h @@ -70,9 +70,9 @@ public: /** - * 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 { @@ -187,14 +187,14 @@ public: * * @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 &target, const double &radius, - std::vector > &matches, - bool sorted=false) const; + std::vector > get_points_within_ball(const Point &target, + const double &radius, + bool sorted=false) const; /** * Fill two vectors with the indices and distances of the closest @@ -204,12 +204,12 @@ public: * 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 &target, - std::vector &indices, - std::vector &distances) const; + std::vector > get_closest_points(const Point &target, + const unsigned int n_points) const; private: /** diff --git a/source/numerics/kdtree_distance.cc b/source/numerics/kdtree_distance.cc index 616a3e847f..72cb1c863b 100644 --- a/source/numerics/kdtree_distance.cc +++ b/source/numerics/kdtree_distance.cc @@ -18,10 +18,11 @@ KDTreeDistance::KDTreeDistance(const unsigned int &max_leaf_size, template -unsigned int KDTreeDistance::get_points_within_ball(const Point ¢er, const double &radius, - std::vector > &matches, - bool sorted) const +std::vector > KDTreeDistance::get_points_within_ball(const Point ¢er, + const double &radius, + bool sorted) const { + std::vector > matches; Assert(adaptor, ExcNotInitialized()); Assert(kdtree, ExcInternalError()); @@ -30,19 +31,24 @@ unsigned int KDTreeDistance::get_points_within_ball(const Point ¢e nanoflann::SearchParams params; params.sorted = sorted; - return kdtree->radiusSearch(¢er[0], radius, matches, params); + kdtree->radiusSearch(¢er[0], radius, matches, params); + return matches; } template -void KDTreeDistance::get_closest_points(const Point &target, - std::vector &indices, - std::vector &distances) const +std::vector > KDTreeDistance::get_closest_points(const Point &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 indices(n_points); + std::vector distances(n_points); + std::vector > matches(n_points); + + kdtree->knnSearch(&target[0], n_points, &indices[0], &distances[0]); + for(unsigned int i=0; i diff --git a/tests/numerics/kdtree_distance_01.cc b/tests/numerics/kdtree_distance_01.cc new file mode 100644 index 0000000000..3946f396dc --- /dev/null +++ b/tests/numerics/kdtree_distance_01.cc @@ -0,0 +1,56 @@ +//----------------------------------------------------------- +// +// 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 + + +int main () +{ + initlog(); + + KDTreeDistance<2> kdtree; + + std::vector > 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 > 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; + } +} diff --git a/tests/numerics/kdtree_distance_01.output b/tests/numerics/kdtree_distance_01.output new file mode 100644 index 0000000000..e1196a586a --- /dev/null +++ b/tests/numerics/kdtree_distance_01.output @@ -0,0 +1,10 @@ + +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 diff --git a/tests/numerics/kdtree_distance_02.cc b/tests/numerics/kdtree_distance_02.cc new file mode 100644 index 0000000000..0706c20825 --- /dev/null +++ b/tests/numerics/kdtree_distance_02.cc @@ -0,0 +1,60 @@ +//----------------------------------------------------------- +// +// 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 + +int main () +{ + initlog(); + + KDTreeDistance<2> kdtree; + + std::vector > 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 > 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 indices; + std::vector distances; + + // Get closest points. Do a few rounds + for (auto &p : test_points) + { + for (unsigned int i=1; i + +int main () +{ + initlog(); + + KDTreeDistance<2> kdtree; + + std::vector > 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 > test_points; + test_points.push_back(Point<2>(2, 0)); + test_points.push_back(Point<2>(2, 2)); + + std::vector 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 > 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; + } + } + } +} diff --git a/tests/numerics/kdtree_distance_03.output b/tests/numerics/kdtree_distance_03.output new file mode 100644 index 0000000000..e4e38ed33b --- /dev/null +++ b/tests/numerics/kdtree_distance_03.output @@ -0,0 +1,24 @@ + +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