Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
kdtree.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_numerics_kdtree_h
17 #define dealii_numerics_kdtree_h
18 
19 #include <deal.II/base/config.h>
20 
21 #ifdef DEAL_II_WITH_NANOFLANN
22 
23 # include <deal.II/base/point.h>
24 
25 # include <nanoflann.hpp>
26 
27 # include <memory>
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
64 template <int dim>
65 class DEAL_II_DEPRECATED KDTree
66 {
67 public:
91  KDTree(const unsigned int max_leaf_size = 10,
92  const std::vector<Point<dim>> &pts = {});
93 
94 
100  {
104  using coord_t = double;
105 
106 
111  const std::vector<Point<dim>> &points;
112 
113 
118  PointCloudAdaptor(const std::vector<Point<dim>> &_points);
119 
120 
124  size_t
125  kdtree_get_point_count() const;
126 
127 
131  coord_t
132  kdtree_distance(const coord_t * p1,
133  const std::size_t idx_p2,
134  const std::size_t size) const;
135 
136 
140  coord_t
141  kdtree_get_pt(const std::size_t idx, const int d) const;
142 
143 
151  template <class BBOX>
152  bool
153  kdtree_get_bbox(BBOX &) const;
154  };
155 
156 
160  using NanoFlannKDTree = typename nanoflann::KDTreeSingleIndexAdaptor<
161  nanoflann::L2_Simple_Adaptor<double, PointCloudAdaptor>,
163  dim,
164  unsigned int>;
165 
166 
184  void
185  set_points(const std::vector<Point<dim>> &pts);
186 
187 
191  const Point<dim> &operator[](const unsigned int i) const;
192 
193 
197  unsigned int
198  size() const;
199 
200 
214  std::vector<std::pair<unsigned int, double>>
215  get_points_within_ball(const Point<dim> &target,
216  const double radius,
217  const bool sorted = false) const;
218 
228  std::vector<std::pair<unsigned int, double>>
229  get_closest_points(const Point<dim> & target,
230  const unsigned int n_points) const;
231 
232 private:
236  const unsigned int max_leaf_size;
237 
238 
242  std::unique_ptr<PointCloudAdaptor> adaptor;
243 
244 
248  std::unique_ptr<NanoFlannKDTree> kdtree;
249 };
250 
251 
252 //------------ inline functions -------------
253 # ifndef DOXYGEN
254 
255 template <int dim>
256 inline unsigned int
257 KDTree<dim>::size() const
258 {
259  if (adaptor)
260  return adaptor->points.size();
261  else
262  return 0;
263 }
264 
265 
266 
267 template <int dim>
268 inline const Point<dim> &KDTree<dim>::operator[](const unsigned int i) const
269 {
270  AssertIndexRange(i, size());
271  return adaptor->points[i];
272 }
273 
274 
275 
276 template <int dim>
278  const std::vector<Point<dim>> &_points)
279  : points(_points)
280 {}
281 
282 
283 
284 template <int dim>
285 inline size_t
287 {
288  return points.size();
289 }
290 
291 
292 
293 template <int dim>
294 inline double
296  int d) const
297 {
298  AssertIndexRange(d, dim);
299  return points[idx][d];
300 }
301 
302 
303 
304 template <int dim>
305 template <class BBOX>
306 inline bool
308 {
309  return false;
310 }
311 
312 
313 
314 template <int dim>
315 inline double
317  const std::size_t idx_p2,
318  const std::size_t size) const
319 {
320  AssertDimension(size, dim);
321  double res = 0.0;
322  for (std::size_t d = 0; d < size; ++d)
323  res += (p1[d] - points[idx_p2][d]) * (p1[d] - points[idx_p2][d]);
324  return std::sqrt(res);
325 }
326 # endif
327 
328 DEAL_II_NAMESPACE_CLOSE
329 
330 #endif // DEAL_II_WITH_NANO_FLANN
331 #endif
const std::vector< Point< dim > > & points
Definition: kdtree.h:111
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1567
Definition: kdtree.h:65
coord_t kdtree_distance(const coord_t *p1, const std::size_t idx_p2, const std::size_t size) const
#define AssertIndexRange(index, range)
Definition: exceptions.h:1637
std::unique_ptr< NanoFlannKDTree > kdtree
Definition: kdtree.h:248
typename nanoflann::KDTreeSingleIndexAdaptor< nanoflann::L2_Simple_Adaptor< double, PointCloudAdaptor >, PointCloudAdaptor, dim, unsigned int > NanoFlannKDTree
Definition: kdtree.h:164
PointCloudAdaptor(const std::vector< Point< dim >> &_points)
size_t kdtree_get_point_count() const
const Point< dim > & operator[](const unsigned int i) const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
unsigned int size() const
coord_t kdtree_get_pt(const std::size_t idx, const int d) const
bool kdtree_get_bbox(BBOX &) const
std::unique_ptr< PointCloudAdaptor > adaptor
Definition: kdtree.h:242
const unsigned int max_leaf_size
Definition: kdtree.h:236