ACloudViewer  3.9.4
A Modern Library for 3D Data Processing
NanoFlannImpl.h
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // - CloudViewer: www.cloudViewer.org -
3 // ----------------------------------------------------------------------------
4 // Copyright (c) 2018-2024 www.cloudViewer.org
5 // SPDX-License-Identifier: MIT
6 // ----------------------------------------------------------------------------
7 #pragma once
8 
9 #include <Helper.h>
10 #include <tbb/parallel_for.h>
11 
12 #include <algorithm>
13 #include <mutex>
14 #include <nanoflann.hpp>
15 
19 
20 namespace cloudViewer {
21 namespace core {
22 namespace nns {
23 
25 template <int METRIC, class TReal, class TIndex>
29  struct DataAdaptor {
30  DataAdaptor(size_t dataset_size,
31  int dimension,
32  const TReal *const data_ptr)
33  : dataset_size_(dataset_size),
34  dimension_(dimension),
35  data_ptr_(data_ptr) {}
36 
37  inline size_t kdtree_get_point_count() const { return dataset_size_; }
38 
39  inline TReal kdtree_get_pt(const size_t idx, const size_t dim) const {
40  return data_ptr_[idx * dimension_ + dim];
41  }
42 
43  template <class BBOX>
44  bool kdtree_get_bbox(BBOX &) const {
45  return false;
46  }
47 
48  size_t dataset_size_ = 0;
49  int dimension_ = 0;
50  const TReal *const data_ptr_;
51  };
52 
54  template <int M, typename fake = void>
56 
57  template <typename fake>
58  struct SelectNanoflannAdaptor<L2, fake> {
59  typedef nanoflann::L2_Adaptor<TReal, DataAdaptor, TReal> adaptor_t;
60  };
61 
62  template <typename fake>
63  struct SelectNanoflannAdaptor<L1, fake> {
64  typedef nanoflann::L1_Adaptor<TReal, DataAdaptor, TReal> adaptor_t;
65  };
66 
68  typedef nanoflann::KDTreeSingleIndexAdaptor<
71  -1,
72  TIndex>
74 
75  NanoFlannIndexHolder(size_t dataset_size,
76  int dimension,
77  const TReal *data_ptr) {
78  adaptor_.reset(new DataAdaptor(dataset_size, dimension, data_ptr));
79  index_.reset(new KDTree_t(dimension, *adaptor_.get()));
80  index_->buildIndex();
81  }
82 
83  std::unique_ptr<KDTree_t> index_;
84  std::unique_ptr<DataAdaptor> adaptor_;
85 };
86 namespace impl {
87 
88 namespace {
89 template <class T, class TIndex, int METRIC>
90 void _BuildKdTree(size_t num_points,
91  const T *const points,
92  size_t dimension,
93  NanoFlannIndexHolderBase **holder) {
94  *holder = new NanoFlannIndexHolder<METRIC, T, TIndex>(num_points, dimension,
95  points);
96 }
97 
98 template <class T, class TIndex, class OUTPUT_ALLOCATOR, int METRIC>
99 void _KnnSearchCPU(NanoFlannIndexHolderBase *holder,
100  int64_t *query_neighbors_row_splits,
101  size_t num_points,
102  const T *const points,
103  size_t num_queries,
104  const T *const queries,
105  const size_t dimension,
106  int knn,
107  bool ignore_query_point,
108  bool return_distances,
109  OUTPUT_ALLOCATOR &output_allocator) {
110  // return empty indices array if there are no points
111  if (num_queries == 0 || num_points == 0 || holder == nullptr) {
112  std::fill(query_neighbors_row_splits,
113  query_neighbors_row_splits + num_queries + 1, 0);
114  TIndex *indices_ptr;
115  output_allocator.AllocIndices(&indices_ptr, 0);
116 
117  T *distances_ptr;
118  output_allocator.AllocDistances(&distances_ptr, 0);
119  return;
120  }
121 
122  std::vector<std::vector<TIndex>> neighbors_indices(num_queries);
123  std::vector<std::vector<T>> neighbors_distances(num_queries);
124  std::vector<uint32_t> neighbors_count(num_queries, 0);
125 
126  // cast NanoFlannIndexHolder
127  auto holder_ =
128  static_cast<NanoFlannIndexHolder<METRIC, T, TIndex> *>(holder);
129 
130  tbb::parallel_for(
131  tbb::blocked_range<size_t>(0, num_queries),
132  [&](const tbb::blocked_range<size_t> &r) {
133  std::vector<TIndex> result_indices(knn);
134  std::vector<T> result_distances(knn);
135  for (size_t i = r.begin(); i != r.end(); ++i) {
136  size_t num_valid = holder_->index_->knnSearch(
137  &queries[i * dimension], knn, result_indices.data(),
138  result_distances.data());
139 
140  int num_neighbors = 0;
141  for (size_t valid_i = 0; valid_i < num_valid; ++valid_i) {
142  TIndex idx = result_indices[valid_i];
143  if (ignore_query_point &&
144  std::equal(&queries[i * dimension],
145  &queries[i * dimension] + dimension,
146  &points[idx * dimension])) {
147  continue;
148  }
149  neighbors_indices[i].push_back(idx);
150  if (return_distances) {
151  T dist = result_distances[valid_i];
152  neighbors_distances[i].push_back(dist);
153  }
154  ++num_neighbors;
155  }
156  neighbors_count[i] = num_neighbors;
157  }
158  });
159 
160  query_neighbors_row_splits[0] = 0;
161  utility::InclusivePrefixSum(neighbors_count.data(),
162  neighbors_count.data() + neighbors_count.size(),
163  query_neighbors_row_splits + 1);
164 
165  int64_t num_indices = query_neighbors_row_splits[num_queries];
166 
167  TIndex *indices_ptr;
168  output_allocator.AllocIndices(&indices_ptr, num_indices);
169  T *distances_ptr;
170  if (return_distances)
171  output_allocator.AllocDistances(&distances_ptr, num_indices);
172  else
173  output_allocator.AllocDistances(&distances_ptr, 0);
174 
175  std::fill(neighbors_count.begin(), neighbors_count.end(), 0);
176 
177  // fill output index and distance arrays
178  tbb::parallel_for(tbb::blocked_range<size_t>(0, num_queries),
179  [&](const tbb::blocked_range<size_t> &r) {
180  for (size_t i = r.begin(); i != r.end(); ++i) {
181  int64_t start_idx = query_neighbors_row_splits[i];
182  std::copy(neighbors_indices[i].begin(),
183  neighbors_indices[i].end(),
184  &indices_ptr[start_idx]);
185 
186  if (return_distances) {
187  std::copy(neighbors_distances[i].begin(),
188  neighbors_distances[i].end(),
189  &distances_ptr[start_idx]);
190  }
191  }
192  });
193 }
194 
195 template <class T, class TIndex, class OUTPUT_ALLOCATOR, int METRIC>
196 void _RadiusSearchCPU(NanoFlannIndexHolderBase *holder,
197  int64_t *query_neighbors_row_splits,
198  size_t num_points,
199  const T *const points,
200  size_t num_queries,
201  const T *const queries,
202  const size_t dimension,
203  const T *const radii,
204  bool ignore_query_point,
205  bool return_distances,
206  bool normalize_distances,
207  bool sort,
208  OUTPUT_ALLOCATOR &output_allocator) {
209  if (num_queries == 0 || num_points == 0 || holder == nullptr) {
210  std::fill(query_neighbors_row_splits,
211  query_neighbors_row_splits + num_queries + 1, 0);
212  TIndex *indices_ptr;
213  output_allocator.AllocIndices(&indices_ptr, 0);
214 
215  T *distances_ptr;
216  output_allocator.AllocDistances(&distances_ptr, 0);
217  return;
218  }
219 
220  std::vector<std::vector<TIndex>> neighbors_indices(num_queries);
221  std::vector<std::vector<T>> neighbors_distances(num_queries);
222  std::vector<uint32_t> neighbors_count(num_queries, 0);
223 
224  nanoflann::SearchParameters params;
225  params.sorted = sort;
226 
227  auto holder_ =
228  static_cast<NanoFlannIndexHolder<METRIC, T, TIndex> *>(holder);
229  tbb::parallel_for(
230  tbb::blocked_range<size_t>(0, num_queries),
231  [&](const tbb::blocked_range<size_t> &r) {
232  std::vector<nanoflann::ResultItem<TIndex, T>> search_result;
233  for (size_t i = r.begin(); i != r.end(); ++i) {
234  T radius = radii[i];
235  if (METRIC == L2) {
236  radius = radius * radius;
237  }
238 
239  holder_->index_->radiusSearch(&queries[i * dimension],
240  radius, search_result,
241  params);
242 
243  int num_neighbors = 0;
244  for (const auto &idx_dist : search_result) {
245  if (ignore_query_point &&
246  std::equal(&queries[i * dimension],
247  &queries[i * dimension] + dimension,
248  &points[static_cast<TIndex>(
249  idx_dist.first) *
250  dimension])) {
251  continue;
252  }
253  neighbors_indices[i].push_back(
254  static_cast<TIndex>(idx_dist.first));
255  if (return_distances) {
256  neighbors_distances[i].push_back(idx_dist.second);
257  }
258  ++num_neighbors;
259  }
260  neighbors_count[i] = num_neighbors;
261  }
262  });
263 
264  query_neighbors_row_splits[0] = 0;
265  utility::InclusivePrefixSum(neighbors_count.data(),
266  neighbors_count.data() + neighbors_count.size(),
267  query_neighbors_row_splits + 1);
268 
269  int64_t num_indices = query_neighbors_row_splits[num_queries];
270 
271  TIndex *indices_ptr;
272  output_allocator.AllocIndices(&indices_ptr, num_indices);
273  T *distances_ptr;
274  if (return_distances)
275  output_allocator.AllocDistances(&distances_ptr, num_indices);
276  else
277  output_allocator.AllocDistances(&distances_ptr, 0);
278 
279  std::fill(neighbors_count.begin(), neighbors_count.end(), 0);
280 
281  // fill output index and distance arrays
282  tbb::parallel_for(
283  tbb::blocked_range<size_t>(0, num_queries),
284  [&](const tbb::blocked_range<size_t> &r) {
285  for (size_t i = r.begin(); i != r.end(); ++i) {
286  int64_t start_idx = query_neighbors_row_splits[i];
287  std::copy(neighbors_indices[i].begin(),
288  neighbors_indices[i].end(),
289  &indices_ptr[start_idx]);
290  if (return_distances) {
291  std::transform(neighbors_distances[i].begin(),
292  neighbors_distances[i].end(),
293  &distances_ptr[start_idx], [&](T dist) {
294  T d = dist;
295  if (normalize_distances) {
296  if (METRIC == L2) {
297  d /= (radii[i] * radii[i]);
298  } else {
299  d /= radii[i];
300  }
301  }
302  return d;
303  });
304  }
305  }
306  });
307 }
308 
309 template <class T, class TIndex, class OUTPUT_ALLOCATOR, int METRIC>
310 void _HybridSearchCPU(NanoFlannIndexHolderBase *holder,
311  size_t num_points,
312  const T *const points,
313  size_t num_queries,
314  const T *const queries,
315  const size_t dimension,
316  const T radius,
317  const int max_knn,
318  bool ignore_query_point,
319  bool return_distances,
320  OUTPUT_ALLOCATOR &output_allocator) {
321  if (num_queries == 0 || num_points == 0 || holder == nullptr) {
322  TIndex *indices_ptr, *counts_ptr;
323  output_allocator.AllocIndices(&indices_ptr, 0);
324  output_allocator.AllocCounts(&counts_ptr, 0);
325 
326  T *distances_ptr;
327  output_allocator.AllocDistances(&distances_ptr, 0);
328  return;
329  }
330 
331  T radius_squared = radius * radius;
332  TIndex *indices_ptr, *counts_ptr;
333  T *distances_ptr;
334 
335  size_t num_indices = static_cast<size_t>(max_knn) * num_queries;
336  output_allocator.AllocIndices(&indices_ptr, num_indices);
337  output_allocator.AllocDistances(&distances_ptr, num_indices);
338  output_allocator.AllocCounts(&counts_ptr, num_queries);
339 
340  nanoflann::SearchParameters params;
341  params.sorted = true;
342 
343  auto holder_ =
344  static_cast<NanoFlannIndexHolder<METRIC, T, TIndex> *>(holder);
345  tbb::parallel_for(
346  tbb::blocked_range<size_t>(0, num_queries),
347  [&](const tbb::blocked_range<size_t> &r) {
348  std::vector<nanoflann::ResultItem<TIndex, T>> ret_matches;
349  for (size_t i = r.begin(); i != r.end(); ++i) {
350  size_t num_results = holder_->index_->radiusSearch(
351  &queries[i * dimension], radius_squared,
352  ret_matches, params);
353  ret_matches.resize(num_results);
354 
355  TIndex count_i = static_cast<TIndex>(num_results);
356  count_i = count_i < max_knn ? count_i : max_knn;
357  counts_ptr[i] = count_i;
358 
359  int neighbor_idx = 0;
360  for (auto it = ret_matches.begin();
361  it < ret_matches.end() && neighbor_idx < max_knn;
362  it++, neighbor_idx++) {
363  indices_ptr[i * max_knn + neighbor_idx] =
364  static_cast<TIndex>(it->first);
365  distances_ptr[i * max_knn + neighbor_idx] = it->second;
366  }
367 
368  while (neighbor_idx < max_knn) {
369  indices_ptr[i * max_knn + neighbor_idx] = -1;
370  distances_ptr[i * max_knn + neighbor_idx] = 0;
371  neighbor_idx += 1;
372  }
373  }
374  });
375 }
376 
377 } // namespace
378 
393 template <class T, class TIndex>
394 std::unique_ptr<NanoFlannIndexHolderBase> BuildKdTree(size_t num_points,
395  const T *const points,
396  size_t dimension,
397  const Metric metric) {
398  NanoFlannIndexHolderBase *holder = nullptr;
399 #define FN_PARAMETERS num_points, points, dimension, &holder
400 
401 #define CALL_TEMPLATE(METRIC) \
402  if (METRIC == metric) { \
403  _BuildKdTree<T, TIndex, METRIC>(FN_PARAMETERS); \
404  }
405 
406 #define CALL_TEMPLATE2 \
407  CALL_TEMPLATE(L1) \
408  CALL_TEMPLATE(L2)
409 
411 
412 #undef CALL_TEMPLATE
413 #undef CALL_TEMPLATE2
414 
415 #undef FN_PARAMETERS
416  return std::unique_ptr<NanoFlannIndexHolderBase>(holder);
417 }
418 
471 template <class T, class TIndex, class OUTPUT_ALLOCATOR>
473  int64_t *query_neighbors_row_splits,
474  size_t num_points,
475  const T *const points,
476  size_t num_queries,
477  const T *const queries,
478  const size_t dimension,
479  int knn,
480  const Metric metric,
481  bool ignore_query_point,
482  bool return_distances,
483  OUTPUT_ALLOCATOR &output_allocator) {
484 #define FN_PARAMETERS \
485  holder, query_neighbors_row_splits, num_points, points, num_queries, \
486  queries, dimension, knn, ignore_query_point, return_distances, \
487  output_allocator
488 
489 #define CALL_TEMPLATE(METRIC) \
490  if (METRIC == metric) { \
491  _KnnSearchCPU<T, TIndex, OUTPUT_ALLOCATOR, METRIC>(FN_PARAMETERS); \
492  }
493 
494 #define CALL_TEMPLATE2 \
495  CALL_TEMPLATE(L1) \
496  CALL_TEMPLATE(L2)
497 
499 
500 #undef CALL_TEMPLATE
501 #undef CALL_TEMPLATE2
502 
503 #undef FN_PARAMETERS
504 }
505 
565 template <class T, class TIndex, class OUTPUT_ALLOCATOR>
567  int64_t *query_neighbors_row_splits,
568  size_t num_points,
569  const T *const points,
570  size_t num_queries,
571  const T *const queries,
572  const size_t dimension,
573  const T *const radii,
574  const Metric metric,
575  bool ignore_query_point,
576  bool return_distances,
577  bool normalize_distances,
578  bool sort,
579  OUTPUT_ALLOCATOR &output_allocator) {
580 #define FN_PARAMETERS \
581  holder, query_neighbors_row_splits, num_points, points, num_queries, \
582  queries, dimension, radii, ignore_query_point, return_distances, \
583  normalize_distances, sort, output_allocator
584 
585 #define CALL_TEMPLATE(METRIC) \
586  if (METRIC == metric) { \
587  _RadiusSearchCPU<T, TIndex, OUTPUT_ALLOCATOR, METRIC>(FN_PARAMETERS); \
588  }
589 
590 #define CALL_TEMPLATE2 \
591  CALL_TEMPLATE(L1) \
592  CALL_TEMPLATE(L2)
593 
595 
596 #undef CALL_TEMPLATE
597 #undef CALL_TEMPLATE2
598 
599 #undef FN_PARAMETERS
600 }
601 
653 template <class T, class TIndex, class OUTPUT_ALLOCATOR>
655  size_t num_points,
656  const T *const points,
657  size_t num_queries,
658  const T *const queries,
659  const size_t dimension,
660  const T radius,
661  const int max_knn,
662  const Metric metric,
663  bool ignore_query_point,
664  bool return_distances,
665  OUTPUT_ALLOCATOR &output_allocator) {
666 #define FN_PARAMETERS \
667  holder, num_points, points, num_queries, queries, dimension, radius, \
668  max_knn, ignore_query_point, return_distances, output_allocator
669 
670 #define CALL_TEMPLATE(METRIC) \
671  if (METRIC == metric) { \
672  _HybridSearchCPU<T, TIndex, OUTPUT_ALLOCATOR, METRIC>(FN_PARAMETERS); \
673  }
674 
675 #define CALL_TEMPLATE2 \
676  CALL_TEMPLATE(L1) \
677  CALL_TEMPLATE(L2)
678 
680 
681 #undef CALL_TEMPLATE
682 #undef CALL_TEMPLATE2
683 
684 #undef FN_PARAMETERS
685 }
686 
687 } // namespace impl
688 } // namespace nns
689 } // namespace core
690 } // namespace cloudViewer
int points
#define CALL_TEMPLATE2
cmdLineReadable * params[]
Helper functions for the ml ops.
static double dist(double x1, double y1, double x2, double y2)
Definition: lsd.c:207
@ BBOX
Definition: CVTypes.h:154
void HybridSearchCPU(NanoFlannIndexHolderBase *holder, size_t num_points, const T *const points, size_t num_queries, const T *const queries, const size_t dimension, const T radius, const int max_knn, const Metric metric, bool ignore_query_point, bool return_distances, OUTPUT_ALLOCATOR &output_allocator)
void RadiusSearchCPU(NanoFlannIndexHolderBase *holder, int64_t *query_neighbors_row_splits, size_t num_points, const T *const points, size_t num_queries, const T *const queries, const size_t dimension, const T *const radii, const Metric metric, bool ignore_query_point, bool return_distances, bool normalize_distances, bool sort, OUTPUT_ALLOCATOR &output_allocator)
void KnnSearchCPU(NanoFlannIndexHolderBase *holder, int64_t *query_neighbors_row_splits, size_t num_points, const T *const points, size_t num_queries, const T *const queries, const size_t dimension, int knn, const Metric metric, bool ignore_query_point, bool return_distances, OUTPUT_ALLOCATOR &output_allocator)
std::unique_ptr< NanoFlannIndexHolderBase > BuildKdTree(size_t num_points, const T *const points, size_t dimension, const Metric metric)
void InclusivePrefixSum(const Tin *first, const Tin *last, Tout *out)
Definition: ParallelScan.h:66
Generic file read and write utility for python interface.
std::vector< PointCoordinateType > radii
Definition: qM3C2Tools.cpp:42
Base struct for NanoFlann index holder.
TReal kdtree_get_pt(const size_t idx, const size_t dim) const
Definition: NanoFlannImpl.h:39
DataAdaptor(size_t dataset_size, int dimension, const TReal *const data_ptr)
Definition: NanoFlannImpl.h:30
NanoFlannIndexHolder(size_t dataset_size, int dimension, const TReal *data_ptr)
Definition: NanoFlannImpl.h:75
nanoflann::KDTreeSingleIndexAdaptor< typename SelectNanoflannAdaptor< METRIC >::adaptor_t, DataAdaptor, -1, TIndex > KDTree_t
typedef for KDtree.
Definition: NanoFlannImpl.h:73
std::unique_ptr< DataAdaptor > adaptor_
Definition: NanoFlannImpl.h:84