1 #include <deal.II/numerics/kdtree.h> 3 #ifdef DEAL_II_WITH_NANOFLANN 5 #include<deal.II/base/std_cxx14/memory.h> 14 max_leaf_size(max_leaf_size)
23 std::vector<std::pair<unsigned int, double> >
32 ExcMessage(
"Radius is expected to be positive."));
34 nanoflann::SearchParams params;
35 params.sorted = sorted;
37 std::vector<std::pair<unsigned int, double> > matches;
38 kdtree->radiusSearch(¢er[0], radius, matches, params);
46 std::vector<std::pair<unsigned int, double> >
48 const unsigned int n_points)
const 54 std::vector<unsigned int> indices(n_points);
55 std::vector<double> distances(n_points);
57 kdtree->knnSearch(&target[0], n_points, &indices[0], &distances[0]);
60 std::vector<std::pair<unsigned int, double> > matches(n_points);
61 for (
unsigned int i=0; i<n_points; ++i)
62 matches[i] = std::make_pair(indices[i], distances[i]);
73 adaptor = std_cxx14::make_unique<PointCloudAdaptor>(pts);
74 kdtree = std_cxx14::make_unique<NanoFlannKDTree>(dim,
76 nanoflann::KDTreeSingleIndexAdaptorParams(max_leaf_size));
85 DEAL_II_NAMESPACE_CLOSE
void set_points(const std::vector< Point< dim > > &pts)
static ::ExceptionBase & ExcNotInitialized()
std::vector< std::pair< unsigned int, double > > get_closest_points(const Point< dim > &target, const unsigned int n_points) const
static ::ExceptionBase & ExcMessage(std::string arg1)
KDTree(const unsigned int &max_leaf_size=10, const std::vector< Point< dim > > &pts=std::vector< Point< dim > >())
#define Assert(cond, exc)
std::vector< std::pair< unsigned int, double > > get_points_within_ball(const Point< dim > &target, const double &radius, const bool sorted=false) const
static ::ExceptionBase & ExcInternalError()