nanoflann
C++ header-only ANN library
Loading...
Searching...
No Matches
nanoflann.hpp
1/***********************************************************************
2 * Software License Agreement (BSD License)
3 *
4 * Copyright 2008-2009 Marius Muja (mariusm@cs.ubc.ca). All rights reserved.
5 * Copyright 2008-2009 David G. Lowe (lowe@cs.ubc.ca). All rights reserved.
6 * Copyright 2011-2025 Jose Luis Blanco (joseluisblancoc@gmail.com).
7 * All rights reserved.
8 *
9 * THE BSD LICENSE
10 *
11 * Redistribution and use in source and binary forms, with or without
12 * modification, are permitted provided that the following conditions
13 * are met:
14 *
15 * 1. Redistributions of source code must retain the above copyright
16 * notice, this list of conditions and the following disclaimer.
17 * 2. Redistributions in binary form must reproduce the above copyright
18 * notice, this list of conditions and the following disclaimer in the
19 * documentation and/or other materials provided with the distribution.
20 *
21 * THIS SOFTWARE IS PROVIDED BY THE AUTHOR ``AS IS'' AND ANY EXPRESS OR
22 * IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES
23 * OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
24 * IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR ANY DIRECT, INDIRECT,
25 * INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT
26 * NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
27 * DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
28 * THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
29 * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF
30 * THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
31 *************************************************************************/
32
45#pragma once
46
47#include <algorithm>
48#include <array>
49#include <atomic>
50#include <cassert>
51#include <cmath> // for abs()
52#include <cstdint>
53#include <cstdlib> // for abs()
54#include <functional> // std::reference_wrapper
55#include <future>
56#include <istream>
57#include <limits> // std::numeric_limits
58#include <ostream>
59#include <stack>
60#include <stdexcept>
61#include <unordered_set>
62#include <vector>
63
65#define NANOFLANN_VERSION 0x180
66
67// Avoid conflicting declaration of min/max macros in Windows headers
68#if !defined(NOMINMAX) && \
69 (defined(_WIN32) || defined(_WIN32_) || defined(WIN32) || defined(_WIN64))
70#define NOMINMAX
71#ifdef max
72#undef max
73#undef min
74#endif
75#endif
76// Avoid conflicts with X11 headers
77#ifdef None
78#undef None
79#endif
80
81namespace nanoflann
82{
87template <typename T>
89{
90 return static_cast<T>(3.14159265358979323846);
91}
92
97template <typename T, typename = int>
98struct has_resize : std::false_type
99{
100};
101
102template <typename T>
103struct has_resize<T, decltype((void)std::declval<T>().resize(1), 0)>
104 : std::true_type
105{
106};
107
108template <typename T, typename = int>
109struct has_assign : std::false_type
110{
111};
112
113template <typename T>
114struct has_assign<T, decltype((void)std::declval<T>().assign(1, 0), 0)>
115 : std::true_type
116{
117};
118
122template <typename Container>
123inline typename std::enable_if<has_resize<Container>::value, void>::type resize(
124 Container& c, const size_t nElements)
125{
126 c.resize(nElements);
127}
128
133template <typename Container>
134inline typename std::enable_if<!has_resize<Container>::value, void>::type
135 resize(Container& c, const size_t nElements)
136{
137 if (nElements != c.size())
138 throw std::logic_error("Try to change the size of a std::array.");
139}
140
144template <typename Container, typename T>
145inline typename std::enable_if<has_assign<Container>::value, void>::type assign(
146 Container& c, const size_t nElements, const T& value)
147{
148 c.assign(nElements, value);
149}
150
154template <typename Container, typename T>
155inline typename std::enable_if<!has_assign<Container>::value, void>::type
156 assign(Container& c, const size_t nElements, const T& value)
157{
158 for (size_t i = 0; i < nElements; i++) c[i] = value;
159}
160
163{
165 template <typename PairType>
166 bool operator()(const PairType& p1, const PairType& p2) const
167 {
168 return p1.second < p2.second;
169 }
170};
171
180template <typename IndexType = size_t, typename DistanceType = double>
182{
183 ResultItem() = default;
184 ResultItem(const IndexType index, const DistanceType distance)
185 : first(index), second(distance)
186 {
187 }
188
189 IndexType first;
190 DistanceType second;
191};
192
197template <
198 typename _DistanceType, typename _IndexType = size_t,
199 typename _CountType = size_t>
201{
202 public:
203 using DistanceType = _DistanceType;
204 using IndexType = _IndexType;
205 using CountType = _CountType;
206
207 private:
208 IndexType* indices;
209 DistanceType* dists;
210 CountType capacity;
211 CountType count;
212
213 public:
214 explicit KNNResultSet(CountType capacity_)
215 : indices(nullptr), dists(nullptr), capacity(capacity_), count(0)
216 {
217 }
218
219 void init(IndexType* indices_, DistanceType* dists_)
220 {
221 indices = indices_;
222 dists = dists_;
223 count = 0;
224 }
225
226 CountType size() const { return count; }
227 bool empty() const { return count == 0; }
228 bool full() const { return count == capacity; }
229
235 bool addPoint(DistanceType dist, IndexType index)
236 {
237 CountType i;
238 for (i = count; i > 0; --i)
239 {
242#ifdef NANOFLANN_FIRST_MATCH
243 if ((dists[i - 1] > dist) ||
244 ((dist == dists[i - 1]) && (indices[i - 1] > index)))
245 {
246#else
247 if (dists[i - 1] > dist)
248 {
249#endif
250 if (i < capacity)
251 {
252 dists[i] = dists[i - 1];
253 indices[i] = indices[i - 1];
254 }
255 }
256 else
257 break;
258 }
259 if (i < capacity)
260 {
261 dists[i] = dist;
262 indices[i] = index;
263 }
264 if (count < capacity) count++;
265
266 // tell caller that the search shall continue
267 return true;
268 }
269
272 DistanceType worstDist() const
273 {
274 return (count < capacity || !count)
275 ? std::numeric_limits<DistanceType>::max()
276 : dists[count - 1];
277 }
278
279 void sort()
280 {
281 // already sorted
282 }
283};
284
286template <
287 typename _DistanceType, typename _IndexType = size_t,
288 typename _CountType = size_t>
290{
291 public:
292 using DistanceType = _DistanceType;
293 using IndexType = _IndexType;
294 using CountType = _CountType;
295
296 private:
297 IndexType* indices;
298 DistanceType* dists;
299 CountType capacity;
300 CountType count;
301 DistanceType maximumSearchDistanceSquared;
302
303 public:
304 explicit RKNNResultSet(
305 CountType capacity_, DistanceType maximumSearchDistanceSquared_)
306 : indices(nullptr),
307 dists(nullptr),
308 capacity(capacity_),
309 count(0),
310 maximumSearchDistanceSquared(maximumSearchDistanceSquared_)
311 {
312 }
313
314 void init(IndexType* indices_, DistanceType* dists_)
315 {
316 indices = indices_;
317 dists = dists_;
318 count = 0;
319 if (capacity) dists[capacity - 1] = maximumSearchDistanceSquared;
320 }
321
322 CountType size() const { return count; }
323 bool empty() const { return count == 0; }
324 bool full() const { return count == capacity; }
325
331 bool addPoint(DistanceType dist, IndexType index)
332 {
333 CountType i;
334 for (i = count; i > 0; --i)
335 {
338#ifdef NANOFLANN_FIRST_MATCH
339 if ((dists[i - 1] > dist) ||
340 ((dist == dists[i - 1]) && (indices[i - 1] > index)))
341 {
342#else
343 if (dists[i - 1] > dist)
344 {
345#endif
346 if (i < capacity)
347 {
348 dists[i] = dists[i - 1];
349 indices[i] = indices[i - 1];
350 }
351 }
352 else
353 break;
354 }
355 if (i < capacity)
356 {
357 dists[i] = dist;
358 indices[i] = index;
359 }
360 if (count < capacity) count++;
361
362 // tell caller that the search shall continue
363 return true;
364 }
365
368 DistanceType worstDist() const
369 {
370 return (count < capacity || !count) ? maximumSearchDistanceSquared
371 : dists[count - 1];
372 }
373
374 void sort()
375 {
376 // already sorted
377 }
378};
379
383template <typename _DistanceType, typename _IndexType = size_t>
385{
386 public:
387 using DistanceType = _DistanceType;
388 using IndexType = _IndexType;
389
390 public:
391 const DistanceType radius;
392
393 std::vector<ResultItem<IndexType, DistanceType>>& m_indices_dists;
394
395 explicit RadiusResultSet(
396 DistanceType radius_,
397 std::vector<ResultItem<IndexType, DistanceType>>& indices_dists)
398 : radius(radius_), m_indices_dists(indices_dists)
399 {
400 init();
401 }
402
403 void init() { clear(); }
404 void clear() { m_indices_dists.clear(); }
405
406 size_t size() const { return m_indices_dists.size(); }
407 size_t empty() const { return m_indices_dists.empty(); }
408
409 bool full() const { return true; }
410
416 bool addPoint(DistanceType dist, IndexType index)
417 {
418 if (dist < radius) m_indices_dists.emplace_back(index, dist);
419 return true;
420 }
421
422 DistanceType worstDist() const { return radius; }
423
429 {
430 if (m_indices_dists.empty())
431 throw std::runtime_error(
432 "Cannot invoke RadiusResultSet::worst_item() on "
433 "an empty list of results.");
434 auto it = std::max_element(
435 m_indices_dists.begin(), m_indices_dists.end(), IndexDist_Sorter());
436 return *it;
437 }
438
439 void sort()
440 {
441 std::sort(
442 m_indices_dists.begin(), m_indices_dists.end(), IndexDist_Sorter());
443 }
444};
445
450template <typename T>
451void save_value(std::ostream& stream, const T& value)
452{
453 stream.write(reinterpret_cast<const char*>(&value), sizeof(T));
454}
455
456template <typename T>
457void save_value(std::ostream& stream, const std::vector<T>& value)
458{
459 size_t size = value.size();
460 stream.write(reinterpret_cast<const char*>(&size), sizeof(size_t));
461 stream.write(reinterpret_cast<const char*>(value.data()), sizeof(T) * size);
462}
463
464template <typename T>
465void load_value(std::istream& stream, T& value)
466{
467 stream.read(reinterpret_cast<char*>(&value), sizeof(T));
468}
469
470template <typename T>
471void load_value(std::istream& stream, std::vector<T>& value)
472{
473 size_t size;
474 stream.read(reinterpret_cast<char*>(&size), sizeof(size_t));
475 value.resize(size);
476 stream.read(reinterpret_cast<char*>(value.data()), sizeof(T) * size);
477}
483struct Metric
484{
485};
486
497template <
498 class T, class DataSource, typename _DistanceType = T,
499 typename IndexType = size_t>
501{
502 using ElementType = T;
503 using DistanceType = _DistanceType;
504
505 const DataSource& data_source;
506
507 L1_Adaptor(const DataSource& _data_source) : data_source(_data_source) {}
508
509 DistanceType evalMetric(
510 const T* a, const IndexType b_idx, size_t size,
511 DistanceType worst_dist = -1) const
512 {
513 DistanceType result = DistanceType();
514 const T* last = a + size;
515 const T* lastgroup = last - 3;
516 size_t d = 0;
517
518 /* Process 4 items with each loop for efficiency. */
519 while (a < lastgroup)
520 {
521 const DistanceType diff0 =
522 std::abs(a[0] - data_source.kdtree_get_pt(b_idx, d++));
523 const DistanceType diff1 =
524 std::abs(a[1] - data_source.kdtree_get_pt(b_idx, d++));
525 const DistanceType diff2 =
526 std::abs(a[2] - data_source.kdtree_get_pt(b_idx, d++));
527 const DistanceType diff3 =
528 std::abs(a[3] - data_source.kdtree_get_pt(b_idx, d++));
529 result += diff0 + diff1 + diff2 + diff3;
530 a += 4;
531 if ((worst_dist > 0) && (result > worst_dist)) { return result; }
532 }
533 /* Process last 0-3 components. Not needed for standard vector lengths.
534 */
535 while (a < last)
536 {
537 result += std::abs(*a++ - data_source.kdtree_get_pt(b_idx, d++));
538 }
539 return result;
540 }
541
542 template <typename U, typename V>
543 DistanceType accum_dist(const U a, const V b, const size_t) const
544 {
545 return std::abs(a - b);
546 }
547};
548
559template <
560 class T, class DataSource, typename _DistanceType = T,
561 typename IndexType = size_t>
563{
564 using ElementType = T;
565 using DistanceType = _DistanceType;
566
567 const DataSource& data_source;
568
569 L2_Adaptor(const DataSource& _data_source) : data_source(_data_source) {}
570
571 DistanceType evalMetric(
572 const T* a, const IndexType b_idx, size_t size,
573 DistanceType worst_dist = -1) const
574 {
575 DistanceType result = DistanceType();
576 const T* last = a + size;
577 const T* lastgroup = last - 3;
578 size_t d = 0;
579
580 /* Process 4 items with each loop for efficiency. */
581 while (a < lastgroup)
582 {
583 const DistanceType diff0 =
584 a[0] - data_source.kdtree_get_pt(b_idx, d++);
585 const DistanceType diff1 =
586 a[1] - data_source.kdtree_get_pt(b_idx, d++);
587 const DistanceType diff2 =
588 a[2] - data_source.kdtree_get_pt(b_idx, d++);
589 const DistanceType diff3 =
590 a[3] - data_source.kdtree_get_pt(b_idx, d++);
591 result +=
592 diff0 * diff0 + diff1 * diff1 + diff2 * diff2 + diff3 * diff3;
593 a += 4;
594 if ((worst_dist > 0) && (result > worst_dist)) { return result; }
595 }
596 /* Process last 0-3 components. Not needed for standard vector lengths.
597 */
598 while (a < last)
599 {
600 const DistanceType diff0 =
601 *a++ - data_source.kdtree_get_pt(b_idx, d++);
602 result += diff0 * diff0;
603 }
604 return result;
605 }
606
607 template <typename U, typename V>
608 DistanceType accum_dist(const U a, const V b, const size_t) const
609 {
610 return (a - b) * (a - b);
611 }
612};
613
624template <
625 class T, class DataSource, typename _DistanceType = T,
626 typename IndexType = size_t>
628{
629 using ElementType = T;
630 using DistanceType = _DistanceType;
631
632 const DataSource& data_source;
633
634 L2_Simple_Adaptor(const DataSource& _data_source)
635 : data_source(_data_source)
636 {
637 }
638
639 DistanceType evalMetric(
640 const T* a, const IndexType b_idx, size_t size) const
641 {
642 DistanceType result = DistanceType();
643 for (size_t i = 0; i < size; ++i)
644 {
645 const DistanceType diff =
646 a[i] - data_source.kdtree_get_pt(b_idx, i);
647 result += diff * diff;
648 }
649 return result;
650 }
651
652 template <typename U, typename V>
653 DistanceType accum_dist(const U a, const V b, const size_t) const
654 {
655 return (a - b) * (a - b);
656 }
657};
658
669template <
670 class T, class DataSource, typename _DistanceType = T,
671 typename IndexType = size_t>
673{
674 using ElementType = T;
675 using DistanceType = _DistanceType;
676
677 const DataSource& data_source;
678
679 SO2_Adaptor(const DataSource& _data_source) : data_source(_data_source) {}
680
681 DistanceType evalMetric(
682 const T* a, const IndexType b_idx, size_t size) const
683 {
684 return accum_dist(
685 a[size - 1], data_source.kdtree_get_pt(b_idx, size - 1), size - 1);
686 }
687
690 template <typename U, typename V>
691 DistanceType accum_dist(const U a, const V b, const size_t) const
692 {
693 DistanceType result = DistanceType();
694 DistanceType PI = pi_const<DistanceType>();
695 result = b - a;
696 if (result > PI)
697 result -= 2 * PI;
698 else if (result < -PI)
699 result += 2 * PI;
700 return result;
701 }
702};
703
714template <
715 class T, class DataSource, typename _DistanceType = T,
716 typename IndexType = size_t>
718{
719 using ElementType = T;
720 using DistanceType = _DistanceType;
721
723 distance_L2_Simple;
724
725 SO3_Adaptor(const DataSource& _data_source)
726 : distance_L2_Simple(_data_source)
727 {
728 }
729
730 DistanceType evalMetric(
731 const T* a, const IndexType b_idx, size_t size) const
732 {
733 return distance_L2_Simple.evalMetric(a, b_idx, size);
734 }
735
736 template <typename U, typename V>
737 DistanceType accum_dist(const U a, const V b, const size_t idx) const
738 {
739 return distance_L2_Simple.accum_dist(a, b, idx);
740 }
741};
742
744struct metric_L1 : public Metric
745{
746 template <class T, class DataSource, typename IndexType = size_t>
751};
754struct metric_L2 : public Metric
755{
756 template <class T, class DataSource, typename IndexType = size_t>
761};
765{
766 template <class T, class DataSource, typename IndexType = size_t>
771};
773struct metric_SO2 : public Metric
774{
775 template <class T, class DataSource, typename IndexType = size_t>
780};
782struct metric_SO3 : public Metric
783{
784 template <class T, class DataSource, typename IndexType = size_t>
789};
790
796enum class KDTreeSingleIndexAdaptorFlags
797{
798 None = 0,
799 SkipInitialBuildIndex = 1
800};
801
802inline std::underlying_type<KDTreeSingleIndexAdaptorFlags>::type operator&(
803 KDTreeSingleIndexAdaptorFlags lhs, KDTreeSingleIndexAdaptorFlags rhs)
804{
805 using underlying =
806 typename std::underlying_type<KDTreeSingleIndexAdaptorFlags>::type;
807 return static_cast<underlying>(lhs) & static_cast<underlying>(rhs);
808}
809
812{
814 size_t _leaf_max_size = 10,
815 KDTreeSingleIndexAdaptorFlags _flags =
816 KDTreeSingleIndexAdaptorFlags::None,
817 unsigned int _n_thread_build = 1)
818 : leaf_max_size(_leaf_max_size),
819 flags(_flags),
820 n_thread_build(_n_thread_build)
821 {
822 }
823
824 size_t leaf_max_size;
825 KDTreeSingleIndexAdaptorFlags flags;
826 unsigned int n_thread_build;
827};
828
831{
832 SearchParameters(float eps_ = 0, bool sorted_ = true)
833 : eps(eps_), sorted(sorted_)
834 {
835 }
836
837 float eps;
838 bool sorted;
840};
861{
862 static constexpr size_t WORDSIZE = 16; // WORDSIZE must >= 8
863 static constexpr size_t BLOCKSIZE = 8192;
864
865 /* We maintain memory alignment to word boundaries by requiring that all
866 allocations be in multiples of the machine wordsize. */
867 /* Size of machine word in bytes. Must be power of 2. */
868 /* Minimum number of bytes requested at a time from the system. Must be
869 * multiple of WORDSIZE. */
870
871 using Size = size_t;
872
873 Size remaining_ = 0;
874 void* base_ = nullptr;
875 void* loc_ = nullptr;
876
877 void internal_init()
878 {
879 remaining_ = 0;
880 base_ = nullptr;
881 usedMemory = 0;
882 wastedMemory = 0;
883 }
884
885 public:
886 Size usedMemory = 0;
887 Size wastedMemory = 0;
888
892 PooledAllocator() { internal_init(); }
893
897 ~PooledAllocator() { free_all(); }
898
900 void free_all()
901 {
902 while (base_ != nullptr)
903 {
904 // Get pointer to prev block
905 void* prev = *(static_cast<void**>(base_));
906 ::free(base_);
907 base_ = prev;
908 }
909 internal_init();
910 }
911
916 void* malloc(const size_t req_size)
917 {
918 /* Round size up to a multiple of wordsize. The following expression
919 only works for WORDSIZE that is a power of 2, by masking last bits
920 of incremented size to zero.
921 */
922 const Size size = (req_size + (WORDSIZE - 1)) & ~(WORDSIZE - 1);
923
924 /* Check whether a new block must be allocated. Note that the first
925 word of a block is reserved for a pointer to the previous block.
926 */
927 if (size > remaining_)
928 {
929 wastedMemory += remaining_;
930
931 /* Allocate new storage. */
932 const Size blocksize =
933 size > BLOCKSIZE ? size + WORDSIZE : BLOCKSIZE + WORDSIZE;
934
935 // use the standard C malloc to allocate memory
936 void* m = ::malloc(blocksize);
937 if (!m) { throw std::bad_alloc(); }
938
939 /* Fill first word of new block with pointer to previous block. */
940 static_cast<void**>(m)[0] = base_;
941 base_ = m;
942
943 remaining_ = blocksize - WORDSIZE;
944 loc_ = static_cast<char*>(m) + WORDSIZE;
945 }
946 void* rloc = loc_;
947 loc_ = static_cast<char*>(loc_) + size;
948 remaining_ -= size;
949
950 usedMemory += size;
951
952 return rloc;
953 }
954
962 template <typename T>
963 T* allocate(const size_t count = 1)
964 {
965 T* mem = static_cast<T*>(this->malloc(sizeof(T) * count));
966 return mem;
967 }
968};
977template <int32_t DIM, typename T>
979{
980 using type = std::array<T, DIM>;
981};
983template <typename T>
984struct array_or_vector<-1, T>
985{
986 using type = std::vector<T>;
987};
988
1005template <
1006 class Derived, typename Distance, class DatasetAdaptor, int32_t DIM = -1,
1007 typename index_t = uint32_t>
1009{
1010 public:
1013 void freeIndex(Derived& obj)
1014 {
1015 obj.pool_.free_all();
1016 obj.root_node_ = nullptr;
1017 obj.size_at_index_build_ = 0;
1018 }
1019
1020 using ElementType = typename Distance::ElementType;
1021 using DistanceType = typename Distance::DistanceType;
1022 using IndexType = index_t;
1023
1027 std::vector<IndexType> vAcc_;
1028
1029 using Offset = typename decltype(vAcc_)::size_type;
1030 using Size = typename decltype(vAcc_)::size_type;
1031 using Dimension = int32_t;
1032
1033 /*---------------------------
1034 * Internal Data Structures
1035 * --------------------------*/
1036 struct Node
1037 {
1040 union
1041 {
1042 struct leaf
1043 {
1044 Offset left, right;
1045 } lr;
1046 struct nonleaf
1047 {
1048 Dimension divfeat;
1050 DistanceType divlow, divhigh;
1051 } sub;
1052 } node_type;
1053
1055 Node *child1 = nullptr, *child2 = nullptr;
1056 };
1057
1058 using NodePtr = Node*;
1059 using NodeConstPtr = const Node*;
1060
1062 {
1063 ElementType low, high;
1064 };
1065
1066 NodePtr root_node_ = nullptr;
1067
1068 Size leaf_max_size_ = 0;
1069
1071 Size n_thread_build_ = 1;
1073 Size size_ = 0;
1075 Size size_at_index_build_ = 0;
1076 Dimension dim_ = 0;
1077
1080 using BoundingBox = typename array_or_vector<DIM, Interval>::type;
1081
1084 using distance_vector_t = typename array_or_vector<DIM, DistanceType>::type;
1085
1088
1097
1099 Size size(const Derived& obj) const { return obj.size_; }
1100
1102 Size veclen(const Derived& obj) const { return DIM > 0 ? DIM : obj.dim_; }
1103
1105 ElementType dataset_get(
1106 const Derived& obj, IndexType element, Dimension component) const
1107 {
1108 return obj.dataset_.kdtree_get_pt(element, component);
1109 }
1110
1115 Size usedMemory(const Derived& obj) const
1116 {
1117 return obj.pool_.usedMemory + obj.pool_.wastedMemory +
1118 obj.dataset_.kdtree_get_point_count() *
1119 sizeof(IndexType); // pool memory and vind array memory
1120 }
1121
1126 const Derived& obj, Offset ind, Size count, Dimension element,
1127 ElementType& min_elem, ElementType& max_elem) const
1128 {
1129 min_elem = dataset_get(obj, vAcc_[ind], element);
1130 max_elem = min_elem;
1131 for (Offset i = 1; i < count; ++i)
1132 {
1133 ElementType val = dataset_get(obj, vAcc_[ind + i], element);
1134 if (val < min_elem) min_elem = val;
1135 if (val > max_elem) max_elem = val;
1136 }
1137 }
1138
1149 Derived& obj, const Offset left, const Offset right, BoundingBox& bbox)
1150 {
1151 assert(left < obj.dataset_.kdtree_get_point_count());
1152
1153 NodePtr node = obj.pool_.template allocate<Node>(); // allocate memory
1154 const auto dims = (DIM > 0 ? DIM : obj.dim_);
1155
1156 /* If too few exemplars remain, then make this a leaf node. */
1157 if ((right - left) <= static_cast<Offset>(obj.leaf_max_size_))
1158 {
1159 node->child1 = node->child2 = nullptr; /* Mark as leaf node. */
1160 node->node_type.lr.left = left;
1161 node->node_type.lr.right = right;
1162
1163 // compute bounding-box of leaf points
1164 for (Dimension i = 0; i < dims; ++i)
1165 {
1166 bbox[i].low = dataset_get(obj, obj.vAcc_[left], i);
1167 bbox[i].high = dataset_get(obj, obj.vAcc_[left], i);
1168 }
1169 for (Offset k = left + 1; k < right; ++k)
1170 {
1171 for (Dimension i = 0; i < dims; ++i)
1172 {
1173 const auto val = dataset_get(obj, obj.vAcc_[k], i);
1174 if (bbox[i].low > val) bbox[i].low = val;
1175 if (bbox[i].high < val) bbox[i].high = val;
1176 }
1177 }
1178 }
1179 else
1180 {
1181 /* Determine the index, dimension and value for split plane */
1182 Offset idx;
1183 Dimension cutfeat;
1184 DistanceType cutval;
1185 middleSplit_(obj, left, right - left, idx, cutfeat, cutval, bbox);
1186
1187 node->node_type.sub.divfeat = cutfeat;
1188
1189 /* Recurse on left */
1190 BoundingBox left_bbox(bbox);
1191 left_bbox[cutfeat].high = cutval;
1192 node->child1 = this->divideTree(obj, left, left + idx, left_bbox);
1193
1194 /* Recurse on right */
1195 BoundingBox right_bbox(bbox);
1196 right_bbox[cutfeat].low = cutval;
1197 node->child2 = this->divideTree(obj, left + idx, right, right_bbox);
1198
1199 node->node_type.sub.divlow = left_bbox[cutfeat].high;
1200 node->node_type.sub.divhigh = right_bbox[cutfeat].low;
1201
1202 for (Dimension i = 0; i < dims; ++i)
1203 {
1204 bbox[i].low = std::min(left_bbox[i].low, right_bbox[i].low);
1205 bbox[i].high = std::max(left_bbox[i].high, right_bbox[i].high);
1206 }
1207 }
1208
1209 return node;
1210 }
1211
1225 Derived& obj, const Offset left, const Offset right, BoundingBox& bbox,
1226 std::atomic<unsigned int>& thread_count, std::mutex& mutex)
1227 {
1228 std::unique_lock<std::mutex> lock(mutex);
1229 NodePtr node = obj.pool_.template allocate<Node>(); // allocate memory
1230 lock.unlock();
1231
1232 const auto dims = (DIM > 0 ? DIM : obj.dim_);
1233
1234 /* If too few exemplars remain, then make this a leaf node. */
1235 if ((right - left) <= static_cast<Offset>(obj.leaf_max_size_))
1236 {
1237 node->child1 = node->child2 = nullptr; /* Mark as leaf node. */
1238 node->node_type.lr.left = left;
1239 node->node_type.lr.right = right;
1240
1241 // compute bounding-box of leaf points
1242 for (Dimension i = 0; i < dims; ++i)
1243 {
1244 bbox[i].low = dataset_get(obj, obj.vAcc_[left], i);
1245 bbox[i].high = dataset_get(obj, obj.vAcc_[left], i);
1246 }
1247 for (Offset k = left + 1; k < right; ++k)
1248 {
1249 for (Dimension i = 0; i < dims; ++i)
1250 {
1251 const auto val = dataset_get(obj, obj.vAcc_[k], i);
1252 if (bbox[i].low > val) bbox[i].low = val;
1253 if (bbox[i].high < val) bbox[i].high = val;
1254 }
1255 }
1256 }
1257 else
1258 {
1259 /* Determine the index, dimension and value for split plane */
1260 Offset idx;
1261 Dimension cutfeat;
1262 DistanceType cutval;
1263 middleSplit_(obj, left, right - left, idx, cutfeat, cutval, bbox);
1264
1265 node->node_type.sub.divfeat = cutfeat;
1266
1267 std::future<NodePtr> right_future;
1268
1269 /* Recurse on right concurrently, if possible */
1270
1271 BoundingBox right_bbox(bbox);
1272 right_bbox[cutfeat].low = cutval;
1273 if (++thread_count < n_thread_build_)
1274 {
1275 /* Concurrent thread for right recursion */
1276
1277 right_future = std::async(
1278 std::launch::async, &KDTreeBaseClass::divideTreeConcurrent,
1279 this, std::ref(obj), left + idx, right,
1280 std::ref(right_bbox), std::ref(thread_count),
1281 std::ref(mutex));
1282 }
1283 else { --thread_count; }
1284
1285 /* Recurse on left in this thread */
1286
1287 BoundingBox left_bbox(bbox);
1288 left_bbox[cutfeat].high = cutval;
1289 node->child1 = this->divideTreeConcurrent(
1290 obj, left, left + idx, left_bbox, thread_count, mutex);
1291
1292 if (right_future.valid())
1293 {
1294 /* Block and wait for concurrent right from above */
1295
1296 node->child2 = right_future.get();
1297 --thread_count;
1298 }
1299 else
1300 {
1301 /* Otherwise, recurse on right in this thread */
1302
1303 node->child2 = this->divideTreeConcurrent(
1304 obj, left + idx, right, right_bbox, thread_count, mutex);
1305 }
1306
1307 node->node_type.sub.divlow = left_bbox[cutfeat].high;
1308 node->node_type.sub.divhigh = right_bbox[cutfeat].low;
1309
1310 for (Dimension i = 0; i < dims; ++i)
1311 {
1312 bbox[i].low = std::min(left_bbox[i].low, right_bbox[i].low);
1313 bbox[i].high = std::max(left_bbox[i].high, right_bbox[i].high);
1314 }
1315 }
1316
1317 return node;
1318 }
1319
1320 void middleSplit_(
1321 const Derived& obj, const Offset ind, const Size count, Offset& index,
1322 Dimension& cutfeat, DistanceType& cutval, const BoundingBox& bbox)
1323 {
1324 const auto dims = (DIM > 0 ? DIM : obj.dim_);
1325 const auto EPS = static_cast<DistanceType>(0.00001);
1326 ElementType max_span = bbox[0].high - bbox[0].low;
1327 for (Dimension i = 1; i < dims; ++i)
1328 {
1329 ElementType span = bbox[i].high - bbox[i].low;
1330 if (span > max_span) { max_span = span; }
1331 }
1332 ElementType max_spread = -1;
1333 cutfeat = 0;
1334 ElementType min_elem = 0, max_elem = 0;
1335 for (Dimension i = 0; i < dims; ++i)
1336 {
1337 ElementType span = bbox[i].high - bbox[i].low;
1338 if (span >= (1 - EPS) * max_span)
1339 {
1340 ElementType min_elem_, max_elem_;
1341 computeMinMax(obj, ind, count, i, min_elem_, max_elem_);
1342 ElementType spread = max_elem_ - min_elem_;
1343 if (spread > max_spread)
1344 {
1345 cutfeat = i;
1346 max_spread = spread;
1347 min_elem = min_elem_;
1348 max_elem = max_elem_;
1349 }
1350 }
1351 }
1352 // split in the middle
1353 DistanceType split_val = (bbox[cutfeat].low + bbox[cutfeat].high) / 2;
1354
1355 if (split_val < min_elem)
1356 cutval = min_elem;
1357 else if (split_val > max_elem)
1358 cutval = max_elem;
1359 else
1360 cutval = split_val;
1361
1362 Offset lim1, lim2;
1363 planeSplit(obj, ind, count, cutfeat, cutval, lim1, lim2);
1364
1365 if (lim1 > count / 2)
1366 index = lim1;
1367 else if (lim2 < count / 2)
1368 index = lim2;
1369 else
1370 index = count / 2;
1371 }
1372
1383 const Derived& obj, const Offset ind, const Size count,
1384 const Dimension cutfeat, const DistanceType& cutval, Offset& lim1,
1385 Offset& lim2)
1386 {
1387 /* First pass.
1388 * Determine lim1 with all values less than cutval to the left.
1389 */
1390 Offset left = 0;
1391 Offset right = count - 1;
1392 for (;;)
1393 {
1394 while (left <= right &&
1395 dataset_get(obj, vAcc_[ind + left], cutfeat) < cutval)
1396 ++left;
1397 while (right && left <= right &&
1398 dataset_get(obj, vAcc_[ind + right], cutfeat) >= cutval)
1399 --right;
1400 if (left > right || !right)
1401 break; // "!right" was added to support unsigned Index types
1402 std::swap(vAcc_[ind + left], vAcc_[ind + right]);
1403 ++left;
1404 --right;
1405 }
1406 /* Second pass
1407 * Determine lim2 with all values greater than cutval to the right
1408 * The middle is used for balancing the tree
1409 */
1410 lim1 = left;
1411 right = count - 1;
1412 for (;;)
1413 {
1414 while (left <= right &&
1415 dataset_get(obj, vAcc_[ind + left], cutfeat) <= cutval)
1416 ++left;
1417 while (right && left <= right &&
1418 dataset_get(obj, vAcc_[ind + right], cutfeat) > cutval)
1419 --right;
1420 if (left > right || !right)
1421 break; // "!right" was added to support unsigned Index types
1422 std::swap(vAcc_[ind + left], vAcc_[ind + right]);
1423 ++left;
1424 --right;
1425 }
1426 lim2 = left;
1427 }
1428
1429 DistanceType computeInitialDistances(
1430 const Derived& obj, const ElementType* vec,
1431 distance_vector_t& dists) const
1432 {
1433 assert(vec);
1434 DistanceType dist = DistanceType();
1435
1436 for (Dimension i = 0; i < (DIM > 0 ? DIM : obj.dim_); ++i)
1437 {
1438 if (vec[i] < obj.root_bbox_[i].low)
1439 {
1440 dists[i] =
1441 obj.distance_.accum_dist(vec[i], obj.root_bbox_[i].low, i);
1442 dist += dists[i];
1443 }
1444 if (vec[i] > obj.root_bbox_[i].high)
1445 {
1446 dists[i] =
1447 obj.distance_.accum_dist(vec[i], obj.root_bbox_[i].high, i);
1448 dist += dists[i];
1449 }
1450 }
1451 return dist;
1452 }
1453
1454 static void save_tree(
1455 const Derived& obj, std::ostream& stream, const NodeConstPtr tree)
1456 {
1457 save_value(stream, *tree);
1458 if (tree->child1 != nullptr) { save_tree(obj, stream, tree->child1); }
1459 if (tree->child2 != nullptr) { save_tree(obj, stream, tree->child2); }
1460 }
1461
1462 static void load_tree(Derived& obj, std::istream& stream, NodePtr& tree)
1463 {
1464 tree = obj.pool_.template allocate<Node>();
1465 load_value(stream, *tree);
1466 if (tree->child1 != nullptr) { load_tree(obj, stream, tree->child1); }
1467 if (tree->child2 != nullptr) { load_tree(obj, stream, tree->child2); }
1468 }
1469
1475 void saveIndex(const Derived& obj, std::ostream& stream) const
1476 {
1477 save_value(stream, obj.size_);
1478 save_value(stream, obj.dim_);
1479 save_value(stream, obj.root_bbox_);
1480 save_value(stream, obj.leaf_max_size_);
1481 save_value(stream, obj.vAcc_);
1482 if (obj.root_node_) save_tree(obj, stream, obj.root_node_);
1483 }
1484
1490 void loadIndex(Derived& obj, std::istream& stream)
1491 {
1492 load_value(stream, obj.size_);
1493 load_value(stream, obj.dim_);
1494 load_value(stream, obj.root_bbox_);
1495 load_value(stream, obj.leaf_max_size_);
1496 load_value(stream, obj.vAcc_);
1497 load_tree(obj, stream, obj.root_node_);
1498 }
1499};
1500
1542template <
1543 typename Distance, class DatasetAdaptor, int32_t DIM = -1,
1544 typename index_t = uint32_t>
1546 : public KDTreeBaseClass<
1547 KDTreeSingleIndexAdaptor<Distance, DatasetAdaptor, DIM, index_t>,
1548 Distance, DatasetAdaptor, DIM, index_t>
1549{
1550 public:
1554 Distance, DatasetAdaptor, DIM, index_t>&) = delete;
1555
1557 const DatasetAdaptor& dataset_;
1558
1559 const KDTreeSingleIndexAdaptorParams indexParams;
1560
1561 Distance distance_;
1562
1563 using Base = typename nanoflann::KDTreeBaseClass<
1565 Distance, DatasetAdaptor, DIM, index_t>,
1566 Distance, DatasetAdaptor, DIM, index_t>;
1567
1568 using Offset = typename Base::Offset;
1569 using Size = typename Base::Size;
1570 using Dimension = typename Base::Dimension;
1571
1572 using ElementType = typename Base::ElementType;
1573 using DistanceType = typename Base::DistanceType;
1574 using IndexType = typename Base::IndexType;
1575
1576 using Node = typename Base::Node;
1577 using NodePtr = Node*;
1578
1579 using Interval = typename Base::Interval;
1580
1583 using BoundingBox = typename Base::BoundingBox;
1584
1587 using distance_vector_t = typename Base::distance_vector_t;
1588
1609 template <class... Args>
1611 const Dimension dimensionality, const DatasetAdaptor& inputData,
1612 const KDTreeSingleIndexAdaptorParams& params, Args&&... args)
1613 : dataset_(inputData),
1614 indexParams(params),
1615 distance_(inputData, std::forward<Args>(args)...)
1616 {
1617 init(dimensionality, params);
1618 }
1619
1620 explicit KDTreeSingleIndexAdaptor(
1621 const Dimension dimensionality, const DatasetAdaptor& inputData,
1622 const KDTreeSingleIndexAdaptorParams& params = {})
1623 : dataset_(inputData), indexParams(params), distance_(inputData)
1624 {
1625 init(dimensionality, params);
1626 }
1627
1628 private:
1629 void init(
1630 const Dimension dimensionality,
1631 const KDTreeSingleIndexAdaptorParams& params)
1632 {
1633 Base::size_ = dataset_.kdtree_get_point_count();
1634 Base::size_at_index_build_ = Base::size_;
1635 Base::dim_ = dimensionality;
1636 if (DIM > 0) Base::dim_ = DIM;
1637 Base::leaf_max_size_ = params.leaf_max_size;
1638 if (params.n_thread_build > 0)
1639 {
1640 Base::n_thread_build_ = params.n_thread_build;
1641 }
1642 else
1643 {
1644 Base::n_thread_build_ =
1645 std::max(std::thread::hardware_concurrency(), 1u);
1646 }
1647
1648 if (!(params.flags &
1649 KDTreeSingleIndexAdaptorFlags::SkipInitialBuildIndex))
1650 {
1651 // Build KD-tree:
1652 buildIndex();
1653 }
1654 }
1655
1656 public:
1661 {
1662 Base::size_ = dataset_.kdtree_get_point_count();
1663 Base::size_at_index_build_ = Base::size_;
1664 init_vind();
1665 this->freeIndex(*this);
1666 Base::size_at_index_build_ = Base::size_;
1667 if (Base::size_ == 0) return;
1668 computeBoundingBox(Base::root_bbox_);
1669 // construct the tree
1670 if (Base::n_thread_build_ == 1)
1671 {
1672 Base::root_node_ =
1673 this->divideTree(*this, 0, Base::size_, Base::root_bbox_);
1674 }
1675 else
1676 {
1677#ifndef NANOFLANN_NO_THREADS
1678 std::atomic<unsigned int> thread_count(0u);
1679 std::mutex mutex;
1680 Base::root_node_ = this->divideTreeConcurrent(
1681 *this, 0, Base::size_, Base::root_bbox_, thread_count, mutex);
1682#else /* NANOFLANN_NO_THREADS */
1683 throw std::runtime_error("Multithreading is disabled");
1684#endif /* NANOFLANN_NO_THREADS */
1685 }
1686 }
1687
1707 template <typename RESULTSET>
1709 RESULTSET& result, const ElementType* vec,
1710 const SearchParameters& searchParams = {}) const
1711 {
1712 assert(vec);
1713 if (this->size(*this) == 0) return false;
1714 if (!Base::root_node_)
1715 throw std::runtime_error(
1716 "[nanoflann] findNeighbors() called before building the "
1717 "index.");
1718 float epsError = 1 + searchParams.eps;
1719
1720 // fixed or variable-sized container (depending on DIM)
1721 distance_vector_t dists;
1722 // Fill it with zeros.
1723 auto zero = static_cast<typename RESULTSET::DistanceType>(0);
1724 assign(dists, (DIM > 0 ? DIM : Base::dim_), zero);
1725 DistanceType dist = this->computeInitialDistances(*this, vec, dists);
1726 searchLevel(result, vec, Base::root_node_, dist, dists, epsError);
1727
1728 if (searchParams.sorted) result.sort();
1729
1730 return result.full();
1731 }
1732
1748 template <typename RESULTSET>
1749 Size findWithinBox(RESULTSET& result, const BoundingBox& bbox) const
1750 {
1751 if (this->size(*this) == 0) return 0;
1752 if (!Base::root_node_)
1753 throw std::runtime_error(
1754 "[nanoflann] findWithinBox() called before building the "
1755 "index.");
1756
1757 std::stack<NodePtr> stack;
1758 stack.push(Base::root_node_);
1759
1760 while (!stack.empty())
1761 {
1762 const NodePtr node = stack.top();
1763 stack.pop();
1764
1765 // If this is a leaf node, then do check and return.
1766 if (!node->child1) // (if one node is nullptr, both are)
1767 {
1768 for (Offset i = node->node_type.lr.left;
1769 i < node->node_type.lr.right; ++i)
1770 {
1771 if (contains(bbox, Base::vAcc_[i]))
1772 {
1773 if (!result.addPoint(0, Base::vAcc_[i]))
1774 {
1775 // the resultset doesn't want to receive any more
1776 // points, we're done searching!
1777 return result.size();
1778 }
1779 }
1780 }
1781 }
1782 else
1783 {
1784 const int idx = node->node_type.sub.divfeat;
1785 const auto low_bound = node->node_type.sub.divlow;
1786 const auto high_bound = node->node_type.sub.divhigh;
1787
1788 if (bbox[idx].low <= low_bound) stack.push(node->child1);
1789 if (bbox[idx].high >= high_bound) stack.push(node->child2);
1790 }
1791 }
1792
1793 return result.size();
1794 }
1795
1812 const ElementType* query_point, const Size num_closest,
1813 IndexType* out_indices, DistanceType* out_distances) const
1814 {
1816 resultSet.init(out_indices, out_distances);
1817 findNeighbors(resultSet, query_point);
1818 return resultSet.size();
1819 }
1820
1841 const ElementType* query_point, const DistanceType& radius,
1842 std::vector<ResultItem<IndexType, DistanceType>>& IndicesDists,
1843 const SearchParameters& searchParams = {}) const
1844 {
1846 radius, IndicesDists);
1847 const Size nFound =
1848 radiusSearchCustomCallback(query_point, resultSet, searchParams);
1849 return nFound;
1850 }
1851
1857 template <class SEARCH_CALLBACK>
1859 const ElementType* query_point, SEARCH_CALLBACK& resultSet,
1860 const SearchParameters& searchParams = {}) const
1861 {
1862 findNeighbors(resultSet, query_point, searchParams);
1863 return resultSet.size();
1864 }
1865
1883 const ElementType* query_point, const Size num_closest,
1884 IndexType* out_indices, DistanceType* out_distances,
1885 const DistanceType& radius) const
1886 {
1888 num_closest, radius);
1889 resultSet.init(out_indices, out_distances);
1890 findNeighbors(resultSet, query_point);
1891 return resultSet.size();
1892 }
1893
1896 public:
1900 {
1901 // Create a permutable array of indices to the input vectors.
1902 Base::size_ = dataset_.kdtree_get_point_count();
1903 if (Base::vAcc_.size() != Base::size_) Base::vAcc_.resize(Base::size_);
1904 for (IndexType i = 0; i < static_cast<IndexType>(Base::size_); i++)
1905 Base::vAcc_[i] = i;
1906 }
1907
1908 void computeBoundingBox(BoundingBox& bbox)
1909 {
1910 const auto dims = (DIM > 0 ? DIM : Base::dim_);
1911 resize(bbox, dims);
1912 if (dataset_.kdtree_get_bbox(bbox))
1913 {
1914 // Done! It was implemented in derived class
1915 }
1916 else
1917 {
1918 const Size N = dataset_.kdtree_get_point_count();
1919 if (!N)
1920 throw std::runtime_error(
1921 "[nanoflann] computeBoundingBox() called but "
1922 "no data points found.");
1923 for (Dimension i = 0; i < dims; ++i)
1924 {
1925 bbox[i].low = bbox[i].high =
1926 this->dataset_get(*this, Base::vAcc_[0], i);
1927 }
1928 for (Offset k = 1; k < N; ++k)
1929 {
1930 for (Dimension i = 0; i < dims; ++i)
1931 {
1932 const auto val =
1933 this->dataset_get(*this, Base::vAcc_[k], i);
1934 if (val < bbox[i].low) bbox[i].low = val;
1935 if (val > bbox[i].high) bbox[i].high = val;
1936 }
1937 }
1938 }
1939 }
1940
1941 bool contains(const BoundingBox& bbox, IndexType idx) const
1942 {
1943 const auto dims = (DIM > 0 ? DIM : Base::dim_);
1944 for (Dimension i = 0; i < dims; ++i)
1945 {
1946 const auto point = this->dataset_.kdtree_get_pt(idx, i);
1947 if (point < bbox[i].low || point > bbox[i].high) return false;
1948 }
1949 return true;
1950 }
1951
1958 template <class RESULTSET>
1960 RESULTSET& result_set, const ElementType* vec, const NodePtr node,
1961 DistanceType mindist, distance_vector_t& dists,
1962 const float epsError) const
1963 {
1964 // If this is a leaf node, then do check and return.
1965 if (!node->child1) // (if one node is nullptr, both are)
1966 {
1967 for (Offset i = node->node_type.lr.left;
1968 i < node->node_type.lr.right; ++i)
1969 {
1970 const IndexType accessor = Base::vAcc_[i]; // reorder... : i;
1971 DistanceType dist = distance_.evalMetric(
1972 vec, accessor, (DIM > 0 ? DIM : Base::dim_));
1973 if (dist < result_set.worstDist())
1974 {
1975 if (!result_set.addPoint(dist, Base::vAcc_[i]))
1976 {
1977 // the resultset doesn't want to receive any more
1978 // points, we're done searching!
1979 return false;
1980 }
1981 }
1982 }
1983 return true;
1984 }
1985
1986 /* Which child branch should be taken first? */
1987 Dimension idx = node->node_type.sub.divfeat;
1988 ElementType val = vec[idx];
1989 DistanceType diff1 = val - node->node_type.sub.divlow;
1990 DistanceType diff2 = val - node->node_type.sub.divhigh;
1991
1992 NodePtr bestChild;
1993 NodePtr otherChild;
1994 DistanceType cut_dist;
1995 if ((diff1 + diff2) < 0)
1996 {
1997 bestChild = node->child1;
1998 otherChild = node->child2;
1999 cut_dist =
2000 distance_.accum_dist(val, node->node_type.sub.divhigh, idx);
2001 }
2002 else
2003 {
2004 bestChild = node->child2;
2005 otherChild = node->child1;
2006 cut_dist =
2007 distance_.accum_dist(val, node->node_type.sub.divlow, idx);
2008 }
2009
2010 /* Call recursively to search next level down. */
2011 if (!searchLevel(result_set, vec, bestChild, mindist, dists, epsError))
2012 {
2013 // the resultset doesn't want to receive any more points, we're done
2014 // searching!
2015 return false;
2016 }
2017
2018 DistanceType dst = dists[idx];
2019 mindist = mindist + cut_dist - dst;
2020 dists[idx] = cut_dist;
2021 if (mindist * epsError <= result_set.worstDist())
2022 {
2023 if (!searchLevel(
2024 result_set, vec, otherChild, mindist, dists, epsError))
2025 {
2026 // the resultset doesn't want to receive any more points, we're
2027 // done searching!
2028 return false;
2029 }
2030 }
2031 dists[idx] = dst;
2032 return true;
2033 }
2034
2035 public:
2041 void saveIndex(std::ostream& stream) const
2042 {
2043 Base::saveIndex(*this, stream);
2044 }
2045
2051 void loadIndex(std::istream& stream) { Base::loadIndex(*this, stream); }
2052
2053}; // class KDTree
2054
2092template <
2093 typename Distance, class DatasetAdaptor, int32_t DIM = -1,
2094 typename IndexType = uint32_t>
2096 : public KDTreeBaseClass<
2097 KDTreeSingleIndexDynamicAdaptor_<
2098 Distance, DatasetAdaptor, DIM, IndexType>,
2099 Distance, DatasetAdaptor, DIM, IndexType>
2100{
2101 public:
2105 const DatasetAdaptor& dataset_;
2106
2107 KDTreeSingleIndexAdaptorParams index_params_;
2108
2109 std::vector<int>& treeIndex_;
2110
2111 Distance distance_;
2112
2113 using Base = typename nanoflann::KDTreeBaseClass<
2115 Distance, DatasetAdaptor, DIM, IndexType>,
2116 Distance, DatasetAdaptor, DIM, IndexType>;
2117
2118 using ElementType = typename Base::ElementType;
2119 using DistanceType = typename Base::DistanceType;
2120
2121 using Offset = typename Base::Offset;
2122 using Size = typename Base::Size;
2123 using Dimension = typename Base::Dimension;
2124
2125 using Node = typename Base::Node;
2126 using NodePtr = Node*;
2127
2128 using Interval = typename Base::Interval;
2131 using BoundingBox = typename Base::BoundingBox;
2132
2135 using distance_vector_t = typename Base::distance_vector_t;
2136
2153 const Dimension dimensionality, const DatasetAdaptor& inputData,
2154 std::vector<int>& treeIndex,
2155 const KDTreeSingleIndexAdaptorParams& params =
2157 : dataset_(inputData),
2158 index_params_(params),
2159 treeIndex_(treeIndex),
2160 distance_(inputData)
2161 {
2162 Base::size_ = 0;
2163 Base::size_at_index_build_ = 0;
2164 for (auto& v : Base::root_bbox_) v = {};
2165 Base::dim_ = dimensionality;
2166 if (DIM > 0) Base::dim_ = DIM;
2167 Base::leaf_max_size_ = params.leaf_max_size;
2168 if (params.n_thread_build > 0)
2169 {
2170 Base::n_thread_build_ = params.n_thread_build;
2171 }
2172 else
2173 {
2174 Base::n_thread_build_ =
2175 std::max(std::thread::hardware_concurrency(), 1u);
2176 }
2177 }
2178
2181 const KDTreeSingleIndexDynamicAdaptor_& rhs) = default;
2182
2186 {
2188 std::swap(Base::vAcc_, tmp.Base::vAcc_);
2189 std::swap(Base::leaf_max_size_, tmp.Base::leaf_max_size_);
2190 std::swap(index_params_, tmp.index_params_);
2191 std::swap(treeIndex_, tmp.treeIndex_);
2192 std::swap(Base::size_, tmp.Base::size_);
2193 std::swap(Base::size_at_index_build_, tmp.Base::size_at_index_build_);
2194 std::swap(Base::root_node_, tmp.Base::root_node_);
2195 std::swap(Base::root_bbox_, tmp.Base::root_bbox_);
2196 std::swap(Base::pool_, tmp.Base::pool_);
2197 return *this;
2198 }
2199
2204 {
2205 Base::size_ = Base::vAcc_.size();
2206 this->freeIndex(*this);
2207 Base::size_at_index_build_ = Base::size_;
2208 if (Base::size_ == 0) return;
2209 computeBoundingBox(Base::root_bbox_);
2210 // construct the tree
2211 if (Base::n_thread_build_ == 1)
2212 {
2213 Base::root_node_ =
2214 this->divideTree(*this, 0, Base::size_, Base::root_bbox_);
2215 }
2216 else
2217 {
2218#ifndef NANOFLANN_NO_THREADS
2219 std::atomic<unsigned int> thread_count(0u);
2220 std::mutex mutex;
2221 Base::root_node_ = this->divideTreeConcurrent(
2222 *this, 0, Base::size_, Base::root_bbox_, thread_count, mutex);
2223#else /* NANOFLANN_NO_THREADS */
2224 throw std::runtime_error("Multithreading is disabled");
2225#endif /* NANOFLANN_NO_THREADS */
2226 }
2227 }
2228
2252 template <typename RESULTSET>
2254 RESULTSET& result, const ElementType* vec,
2255 const SearchParameters& searchParams = {}) const
2256 {
2257 assert(vec);
2258 if (this->size(*this) == 0) return false;
2259 if (!Base::root_node_) return false;
2260 float epsError = 1 + searchParams.eps;
2261
2262 // fixed or variable-sized container (depending on DIM)
2263 distance_vector_t dists;
2264 // Fill it with zeros.
2265 assign(
2266 dists, (DIM > 0 ? DIM : Base::dim_),
2267 static_cast<typename distance_vector_t::value_type>(0));
2268 DistanceType dist = this->computeInitialDistances(*this, vec, dists);
2269 searchLevel(result, vec, Base::root_node_, dist, dists, epsError);
2270 return result.full();
2271 }
2272
2288 const ElementType* query_point, const Size num_closest,
2289 IndexType* out_indices, DistanceType* out_distances,
2290 const SearchParameters& searchParams = {}) const
2291 {
2293 resultSet.init(out_indices, out_distances);
2294 findNeighbors(resultSet, query_point, searchParams);
2295 return resultSet.size();
2296 }
2297
2318 const ElementType* query_point, const DistanceType& radius,
2319 std::vector<ResultItem<IndexType, DistanceType>>& IndicesDists,
2320 const SearchParameters& searchParams = {}) const
2321 {
2323 radius, IndicesDists);
2324 const size_t nFound =
2325 radiusSearchCustomCallback(query_point, resultSet, searchParams);
2326 return nFound;
2327 }
2328
2334 template <class SEARCH_CALLBACK>
2336 const ElementType* query_point, SEARCH_CALLBACK& resultSet,
2337 const SearchParameters& searchParams = {}) const
2338 {
2339 findNeighbors(resultSet, query_point, searchParams);
2340 return resultSet.size();
2341 }
2342
2345 public:
2346 void computeBoundingBox(BoundingBox& bbox)
2347 {
2348 const auto dims = (DIM > 0 ? DIM : Base::dim_);
2349 resize(bbox, dims);
2350
2351 if (dataset_.kdtree_get_bbox(bbox))
2352 {
2353 // Done! It was implemented in derived class
2354 }
2355 else
2356 {
2357 const Size N = Base::size_;
2358 if (!N)
2359 throw std::runtime_error(
2360 "[nanoflann] computeBoundingBox() called but "
2361 "no data points found.");
2362 for (Dimension i = 0; i < dims; ++i)
2363 {
2364 bbox[i].low = bbox[i].high =
2365 this->dataset_get(*this, Base::vAcc_[0], i);
2366 }
2367 for (Offset k = 1; k < N; ++k)
2368 {
2369 for (Dimension i = 0; i < dims; ++i)
2370 {
2371 const auto val =
2372 this->dataset_get(*this, Base::vAcc_[k], i);
2373 if (val < bbox[i].low) bbox[i].low = val;
2374 if (val > bbox[i].high) bbox[i].high = val;
2375 }
2376 }
2377 }
2378 }
2379
2384 template <class RESULTSET>
2386 RESULTSET& result_set, const ElementType* vec, const NodePtr node,
2387 DistanceType mindist, distance_vector_t& dists,
2388 const float epsError) const
2389 {
2390 // If this is a leaf node, then do check and return.
2391 if (!node->child1) // (if one node is nullptr, both are)
2392 {
2393 for (Offset i = node->node_type.lr.left;
2394 i < node->node_type.lr.right; ++i)
2395 {
2396 const IndexType index = Base::vAcc_[i]; // reorder... : i;
2397 if (treeIndex_[index] == -1) continue;
2398 DistanceType dist = distance_.evalMetric(
2399 vec, index, (DIM > 0 ? DIM : Base::dim_));
2400 if (dist < result_set.worstDist())
2401 {
2402 if (!result_set.addPoint(
2403 static_cast<typename RESULTSET::DistanceType>(dist),
2404 static_cast<typename RESULTSET::IndexType>(
2405 Base::vAcc_[i])))
2406 {
2407 // the resultset doesn't want to receive any more
2408 // points, we're done searching!
2409 return; // false;
2410 }
2411 }
2412 }
2413 return;
2414 }
2415
2416 /* Which child branch should be taken first? */
2417 Dimension idx = node->node_type.sub.divfeat;
2418 ElementType val = vec[idx];
2419 DistanceType diff1 = val - node->node_type.sub.divlow;
2420 DistanceType diff2 = val - node->node_type.sub.divhigh;
2421
2422 NodePtr bestChild;
2423 NodePtr otherChild;
2424 DistanceType cut_dist;
2425 if ((diff1 + diff2) < 0)
2426 {
2427 bestChild = node->child1;
2428 otherChild = node->child2;
2429 cut_dist =
2430 distance_.accum_dist(val, node->node_type.sub.divhigh, idx);
2431 }
2432 else
2433 {
2434 bestChild = node->child2;
2435 otherChild = node->child1;
2436 cut_dist =
2437 distance_.accum_dist(val, node->node_type.sub.divlow, idx);
2438 }
2439
2440 /* Call recursively to search next level down. */
2441 searchLevel(result_set, vec, bestChild, mindist, dists, epsError);
2442
2443 DistanceType dst = dists[idx];
2444 mindist = mindist + cut_dist - dst;
2445 dists[idx] = cut_dist;
2446 if (mindist * epsError <= result_set.worstDist())
2447 {
2448 searchLevel(result_set, vec, otherChild, mindist, dists, epsError);
2449 }
2450 dists[idx] = dst;
2451 }
2452
2453 public:
2459 void saveIndex(std::ostream& stream) { saveIndex(*this, stream); }
2460
2466 void loadIndex(std::istream& stream) { loadIndex(*this, stream); }
2467};
2468
2483template <
2484 typename Distance, class DatasetAdaptor, int32_t DIM = -1,
2485 typename IndexType = uint32_t>
2487{
2488 public:
2489 using ElementType = typename Distance::ElementType;
2490 using DistanceType = typename Distance::DistanceType;
2491
2492 using Offset = typename KDTreeSingleIndexDynamicAdaptor_<
2493 Distance, DatasetAdaptor, DIM>::Offset;
2494 using Size = typename KDTreeSingleIndexDynamicAdaptor_<
2495 Distance, DatasetAdaptor, DIM>::Size;
2496 using Dimension = typename KDTreeSingleIndexDynamicAdaptor_<
2497 Distance, DatasetAdaptor, DIM>::Dimension;
2498
2499 protected:
2500 Size leaf_max_size_;
2501 Size treeCount_;
2502 Size pointCount_;
2503
2507 const DatasetAdaptor& dataset_;
2508
2511 std::vector<int> treeIndex_;
2512 std::unordered_set<int> removedPoints_;
2513
2514 KDTreeSingleIndexAdaptorParams index_params_;
2515
2516 Dimension dim_;
2517
2519 Distance, DatasetAdaptor, DIM, IndexType>;
2520 std::vector<index_container_t> index_;
2521
2522 public:
2525 const std::vector<index_container_t>& getAllIndices() const
2526 {
2527 return index_;
2528 }
2529
2530 private:
2532 int First0Bit(IndexType num)
2533 {
2534 int pos = 0;
2535 while (num & 1)
2536 {
2537 num = num >> 1;
2538 pos++;
2539 }
2540 return pos;
2541 }
2542
2544 void init()
2545 {
2546 using my_kd_tree_t = KDTreeSingleIndexDynamicAdaptor_<
2547 Distance, DatasetAdaptor, DIM, IndexType>;
2548 std::vector<my_kd_tree_t> index(
2549 treeCount_,
2550 my_kd_tree_t(dim_ /*dim*/, dataset_, treeIndex_, index_params_));
2551 index_ = index;
2552 }
2553
2554 public:
2555 Distance distance_;
2556
2573 const int dimensionality, const DatasetAdaptor& inputData,
2574 const KDTreeSingleIndexAdaptorParams& params =
2576 const size_t maximumPointCount = 1000000000U)
2577 : dataset_(inputData), index_params_(params), distance_(inputData)
2578 {
2579 treeCount_ = static_cast<size_t>(std::log2(maximumPointCount)) + 1;
2580 pointCount_ = 0U;
2581 dim_ = dimensionality;
2582 treeIndex_.clear();
2583 if (DIM > 0) dim_ = DIM;
2584 leaf_max_size_ = params.leaf_max_size;
2585 init();
2586 const size_t num_initial_points = dataset_.kdtree_get_point_count();
2587 if (num_initial_points > 0) { addPoints(0, num_initial_points - 1); }
2588 }
2589
2593 Distance, DatasetAdaptor, DIM, IndexType>&) = delete;
2594
2596 void addPoints(IndexType start, IndexType end)
2597 {
2598 const Size count = end - start + 1;
2599 int maxIndex = 0;
2600 treeIndex_.resize(treeIndex_.size() + count);
2601 for (IndexType idx = start; idx <= end; idx++)
2602 {
2603 const int pos = First0Bit(pointCount_);
2604 maxIndex = std::max(pos, maxIndex);
2605 treeIndex_[pointCount_] = pos;
2606
2607 const auto it = removedPoints_.find(idx);
2608 if (it != removedPoints_.end())
2609 {
2610 removedPoints_.erase(it);
2611 treeIndex_[idx] = pos;
2612 }
2613
2614 for (int i = 0; i < pos; i++)
2615 {
2616 for (int j = 0; j < static_cast<int>(index_[i].vAcc_.size());
2617 j++)
2618 {
2619 index_[pos].vAcc_.push_back(index_[i].vAcc_[j]);
2620 if (treeIndex_[index_[i].vAcc_[j]] != -1)
2621 treeIndex_[index_[i].vAcc_[j]] = pos;
2622 }
2623 index_[i].vAcc_.clear();
2624 }
2625 index_[pos].vAcc_.push_back(idx);
2626 pointCount_++;
2627 }
2628
2629 for (int i = 0; i <= maxIndex; ++i)
2630 {
2631 index_[i].freeIndex(index_[i]);
2632 if (!index_[i].vAcc_.empty()) index_[i].buildIndex();
2633 }
2634 }
2635
2637 void removePoint(size_t idx)
2638 {
2639 if (idx >= pointCount_) return;
2640 removedPoints_.insert(idx);
2641 treeIndex_[idx] = -1;
2642 }
2643
2660 template <typename RESULTSET>
2662 RESULTSET& result, const ElementType* vec,
2663 const SearchParameters& searchParams = {}) const
2664 {
2665 for (size_t i = 0; i < treeCount_; i++)
2666 {
2667 index_[i].findNeighbors(result, &vec[0], searchParams);
2668 }
2669 return result.full();
2670 }
2671};
2672
2698template <
2699 class MatrixType, int32_t DIM = -1, class Distance = nanoflann::metric_L2,
2700 bool row_major = true>
2702{
2703 using self_t =
2705 using num_t = typename MatrixType::Scalar;
2706 using IndexType = typename MatrixType::Index;
2707 using metric_t = typename Distance::template traits<
2708 num_t, self_t, IndexType>::distance_t;
2709
2711 metric_t, self_t,
2712 row_major ? MatrixType::ColsAtCompileTime
2713 : MatrixType::RowsAtCompileTime,
2714 IndexType>;
2715
2716 index_t* index_;
2718
2719 using Offset = typename index_t::Offset;
2720 using Size = typename index_t::Size;
2721 using Dimension = typename index_t::Dimension;
2722
2725 const Dimension dimensionality,
2726 const std::reference_wrapper<const MatrixType>& mat,
2727 const int leaf_max_size = 10, const unsigned int n_thread_build = 1)
2728 : m_data_matrix(mat)
2729 {
2730 const auto dims = row_major ? mat.get().cols() : mat.get().rows();
2731 if (static_cast<Dimension>(dims) != dimensionality)
2732 throw std::runtime_error(
2733 "Error: 'dimensionality' must match column count in data "
2734 "matrix");
2735 if (DIM > 0 && static_cast<int32_t>(dims) != DIM)
2736 throw std::runtime_error(
2737 "Data set dimensionality does not match the 'DIM' template "
2738 "argument");
2739 index_ = new index_t(
2740 dims, *this /* adaptor */,
2742 leaf_max_size, nanoflann::KDTreeSingleIndexAdaptorFlags::None,
2743 n_thread_build));
2744 }
2745
2746 public:
2749
2750 ~KDTreeEigenMatrixAdaptor() { delete index_; }
2751
2752 const std::reference_wrapper<const MatrixType> m_data_matrix;
2753
2762 void query(
2763 const num_t* query_point, const Size num_closest,
2764 IndexType* out_indices, num_t* out_distances) const
2765 {
2766 nanoflann::KNNResultSet<num_t, IndexType> resultSet(num_closest);
2767 resultSet.init(out_indices, out_distances);
2768 index_->findNeighbors(resultSet, query_point);
2769 }
2770
2774 const self_t& derived() const { return *this; }
2775 self_t& derived() { return *this; }
2776
2777 // Must return the number of data points
2778 Size kdtree_get_point_count() const
2779 {
2780 if (row_major)
2781 return m_data_matrix.get().rows();
2782 else
2783 return m_data_matrix.get().cols();
2784 }
2785
2786 // Returns the dim'th component of the idx'th point in the class:
2787 num_t kdtree_get_pt(const IndexType idx, size_t dim) const
2788 {
2789 if (row_major)
2790 return m_data_matrix.get().coeff(idx, IndexType(dim));
2791 else
2792 return m_data_matrix.get().coeff(IndexType(dim), idx);
2793 }
2794
2795 // Optional bounding-box computation: return false to default to a standard
2796 // bbox computation loop.
2797 // Return true if the BBOX was already computed by the class and returned
2798 // in "bb" so it can be avoided to redo it again. Look at bb.size() to
2799 // find out the expected dimensionality (e.g. 2 or 3 for point clouds)
2800 template <class BBOX>
2801 bool kdtree_get_bbox(BBOX& /*bb*/) const
2802 {
2803 return false;
2804 }
2805
2808}; // end of KDTreeEigenMatrixAdaptor
// end of grouping
2812} // namespace nanoflann
Definition nanoflann.hpp:1009
void freeIndex(Derived &obj)
Definition nanoflann.hpp:1013
Size veclen(const Derived &obj) const
Definition nanoflann.hpp:1102
void computeMinMax(const Derived &obj, Offset ind, Size count, Dimension element, ElementType &min_elem, ElementType &max_elem) const
Definition nanoflann.hpp:1125
BoundingBox root_bbox_
Definition nanoflann.hpp:1087
void saveIndex(const Derived &obj, std::ostream &stream) const
Definition nanoflann.hpp:1475
typename array_or_vector< DIM, DistanceType >::type distance_vector_t
Definition nanoflann.hpp:1084
void planeSplit(const Derived &obj, const Offset ind, const Size count, const Dimension cutfeat, const DistanceType &cutval, Offset &lim1, Offset &lim2)
Definition nanoflann.hpp:1382
NodePtr divideTree(Derived &obj, const Offset left, const Offset right, BoundingBox &bbox)
Definition nanoflann.hpp:1148
std::vector< IndexType > vAcc_
Definition nanoflann.hpp:1027
Size size(const Derived &obj) const
Definition nanoflann.hpp:1099
NodePtr divideTreeConcurrent(Derived &obj, const Offset left, const Offset right, BoundingBox &bbox, std::atomic< unsigned int > &thread_count, std::mutex &mutex)
Definition nanoflann.hpp:1224
Size usedMemory(const Derived &obj) const
Definition nanoflann.hpp:1115
void loadIndex(Derived &obj, std::istream &stream)
Definition nanoflann.hpp:1490
PooledAllocator pool_
Definition nanoflann.hpp:1096
ElementType dataset_get(const Derived &obj, IndexType element, Dimension component) const
Helper accessor to the dataset points:
Definition nanoflann.hpp:1105
typename array_or_vector< DIM, Interval >::type BoundingBox
Definition nanoflann.hpp:1080
Definition nanoflann.hpp:1549
bool searchLevel(RESULTSET &result_set, const ElementType *vec, const NodePtr node, DistanceType mindist, distance_vector_t &dists, const float epsError) const
Definition nanoflann.hpp:1959
void saveIndex(std::ostream &stream) const
Definition nanoflann.hpp:2041
Size findWithinBox(RESULTSET &result, const BoundingBox &bbox) const
Definition nanoflann.hpp:1749
Size radiusSearch(const ElementType *query_point, const DistanceType &radius, std::vector< ResultItem< IndexType, DistanceType > > &IndicesDists, const SearchParameters &searchParams={}) const
Definition nanoflann.hpp:1840
void init_vind()
Definition nanoflann.hpp:1899
void buildIndex()
Definition nanoflann.hpp:1660
const DatasetAdaptor & dataset_
Definition nanoflann.hpp:1557
KDTreeSingleIndexAdaptor(const KDTreeSingleIndexAdaptor< Distance, DatasetAdaptor, DIM, index_t > &)=delete
bool findNeighbors(RESULTSET &result, const ElementType *vec, const SearchParameters &searchParams={}) const
Definition nanoflann.hpp:1708
Size rknnSearch(const ElementType *query_point, const Size num_closest, IndexType *out_indices, DistanceType *out_distances, const DistanceType &radius) const
Definition nanoflann.hpp:1882
Size radiusSearchCustomCallback(const ElementType *query_point, SEARCH_CALLBACK &resultSet, const SearchParameters &searchParams={}) const
Definition nanoflann.hpp:1858
typename Base::distance_vector_t distance_vector_t
Definition nanoflann.hpp:1587
void loadIndex(std::istream &stream)
Definition nanoflann.hpp:2051
typename Base::BoundingBox BoundingBox
Definition nanoflann.hpp:1583
KDTreeSingleIndexAdaptor(const Dimension dimensionality, const DatasetAdaptor &inputData, const KDTreeSingleIndexAdaptorParams &params, Args &&... args)
Definition nanoflann.hpp:1610
Size knnSearch(const ElementType *query_point, const Size num_closest, IndexType *out_indices, DistanceType *out_distances) const
Definition nanoflann.hpp:1811
Definition nanoflann.hpp:2100
Size radiusSearch(const ElementType *query_point, const DistanceType &radius, std::vector< ResultItem< IndexType, DistanceType > > &IndicesDists, const SearchParameters &searchParams={}) const
Definition nanoflann.hpp:2317
KDTreeSingleIndexDynamicAdaptor_(const Dimension dimensionality, const DatasetAdaptor &inputData, std::vector< int > &treeIndex, const KDTreeSingleIndexAdaptorParams &params=KDTreeSingleIndexAdaptorParams())
Definition nanoflann.hpp:2152
typename Base::BoundingBox BoundingBox
Definition nanoflann.hpp:2131
const DatasetAdaptor & dataset_
The source of our data.
Definition nanoflann.hpp:2105
Size radiusSearchCustomCallback(const ElementType *query_point, SEARCH_CALLBACK &resultSet, const SearchParameters &searchParams={}) const
Definition nanoflann.hpp:2335
KDTreeSingleIndexDynamicAdaptor_(const KDTreeSingleIndexDynamicAdaptor_ &rhs)=default
void buildIndex()
Definition nanoflann.hpp:2203
void saveIndex(std::ostream &stream)
Definition nanoflann.hpp:2459
typename Base::distance_vector_t distance_vector_t
Definition nanoflann.hpp:2135
void loadIndex(std::istream &stream)
Definition nanoflann.hpp:2466
Size knnSearch(const ElementType *query_point, const Size num_closest, IndexType *out_indices, DistanceType *out_distances, const SearchParameters &searchParams={}) const
Definition nanoflann.hpp:2287
void searchLevel(RESULTSET &result_set, const ElementType *vec, const NodePtr node, DistanceType mindist, distance_vector_t &dists, const float epsError) const
Definition nanoflann.hpp:2385
bool findNeighbors(RESULTSET &result, const ElementType *vec, const SearchParameters &searchParams={}) const
Definition nanoflann.hpp:2253
KDTreeSingleIndexDynamicAdaptor_ operator=(const KDTreeSingleIndexDynamicAdaptor_ &rhs)
Definition nanoflann.hpp:2184
Definition nanoflann.hpp:2487
bool findNeighbors(RESULTSET &result, const ElementType *vec, const SearchParameters &searchParams={}) const
Definition nanoflann.hpp:2661
const DatasetAdaptor & dataset_
The source of our data.
Definition nanoflann.hpp:2507
void removePoint(size_t idx)
Definition nanoflann.hpp:2637
void addPoints(IndexType start, IndexType end)
Definition nanoflann.hpp:2596
KDTreeSingleIndexDynamicAdaptor(const int dimensionality, const DatasetAdaptor &inputData, const KDTreeSingleIndexAdaptorParams &params=KDTreeSingleIndexAdaptorParams(), const size_t maximumPointCount=1000000000U)
Definition nanoflann.hpp:2572
std::vector< int > treeIndex_
Definition nanoflann.hpp:2511
const std::vector< index_container_t > & getAllIndices() const
Definition nanoflann.hpp:2525
Dimension dim_
Dimensionality of each data point.
Definition nanoflann.hpp:2516
KDTreeSingleIndexDynamicAdaptor(const KDTreeSingleIndexDynamicAdaptor< Distance, DatasetAdaptor, DIM, IndexType > &)=delete
Definition nanoflann.hpp:201
bool addPoint(DistanceType dist, IndexType index)
Definition nanoflann.hpp:235
DistanceType worstDist() const
Definition nanoflann.hpp:272
Definition nanoflann.hpp:861
~PooledAllocator()
Definition nanoflann.hpp:897
void free_all()
Definition nanoflann.hpp:900
void * malloc(const size_t req_size)
Definition nanoflann.hpp:916
T * allocate(const size_t count=1)
Definition nanoflann.hpp:963
PooledAllocator()
Definition nanoflann.hpp:892
Definition nanoflann.hpp:290
bool addPoint(DistanceType dist, IndexType index)
Definition nanoflann.hpp:331
DistanceType worstDist() const
Definition nanoflann.hpp:368
Definition nanoflann.hpp:385
ResultItem< IndexType, DistanceType > worst_item() const
Definition nanoflann.hpp:428
bool addPoint(DistanceType dist, IndexType index)
Definition nanoflann.hpp:416
std::enable_if< has_assign< Container >::value, void >::type assign(Container &c, const size_t nElements, const T &value)
Definition nanoflann.hpp:145
T pi_const()
Definition nanoflann.hpp:88
std::enable_if< has_resize< Container >::value, void >::type resize(Container &c, const size_t nElements)
Definition nanoflann.hpp:123
Definition nanoflann.hpp:163
bool operator()(const PairType &p1, const PairType &p2) const
Definition nanoflann.hpp:166
Definition nanoflann.hpp:1062
Definition nanoflann.hpp:1037
DistanceType divlow
The values used for subdivision.
Definition nanoflann.hpp:1050
Offset right
Indices of points in leaf node.
Definition nanoflann.hpp:1044
union nanoflann::KDTreeBaseClass::Node::@0 node_type
Dimension divfeat
Definition nanoflann.hpp:1048
Node * child1
Definition nanoflann.hpp:1055
Definition nanoflann.hpp:2702
void query(const num_t *query_point, const Size num_closest, IndexType *out_indices, num_t *out_distances) const
Definition nanoflann.hpp:2762
KDTreeEigenMatrixAdaptor(const self_t &)=delete
typename index_t::Offset Offset
Definition nanoflann.hpp:2719
KDTreeEigenMatrixAdaptor(const Dimension dimensionality, const std::reference_wrapper< const MatrixType > &mat, const int leaf_max_size=10, const unsigned int n_thread_build=1)
Constructor: takes a const ref to the matrix object with the data points.
Definition nanoflann.hpp:2724
Definition nanoflann.hpp:812
Definition nanoflann.hpp:501
Definition nanoflann.hpp:563
Definition nanoflann.hpp:628
Definition nanoflann.hpp:484
Definition nanoflann.hpp:182
DistanceType second
Distance from sample to query point.
Definition nanoflann.hpp:190
IndexType first
Index of the sample in the dataset.
Definition nanoflann.hpp:189
Definition nanoflann.hpp:673
DistanceType accum_dist(const U a, const V b, const size_t) const
Definition nanoflann.hpp:691
Definition nanoflann.hpp:718
Definition nanoflann.hpp:831
bool sorted
distance (default: true)
Definition nanoflann.hpp:838
float eps
search for eps-approximate neighbours (default: 0)
Definition nanoflann.hpp:837
Definition nanoflann.hpp:979
Definition nanoflann.hpp:110
Definition nanoflann.hpp:99
Definition nanoflann.hpp:748
Definition nanoflann.hpp:745
Definition nanoflann.hpp:758
Definition nanoflann.hpp:768
Definition nanoflann.hpp:765
Definition nanoflann.hpp:755
Definition nanoflann.hpp:777
Definition nanoflann.hpp:774
Definition nanoflann.hpp:786
Definition nanoflann.hpp:783