Go to the documentation of this file.
46 #ifndef NANOFLANN_HPP_
47 #define NANOFLANN_HPP_
59 #if !defined(NOMINMAX) && (defined(_WIN32) || defined(_WIN32_) || defined(WIN32) || defined(_WIN64))
73 #define NANOFLANN_VERSION 0x122
77 template <
typename DistanceType,
typename IndexType =
size_t ,
typename CountType =
size_t >
90 inline void init (IndexType* indices_, DistanceType* dists_)
96 dists [
capacity -1] = (std::numeric_limits<DistanceType>::max)();
99 inline CountType
size ()
const
113 for (i=
count ; i>0; --i) {
114 #ifdef NANOFLANN_FIRST_MATCH // If defined and two points have the same distance, the one with the lowest-index will be returned first.
117 if (
dists [i-1]>dist) {
143 template <
typename DistanceType,
typename IndexType =
size_t >
151 inline RadiusResultSet (DistanceType radius_, std::vector<std::pair<IndexType,DistanceType> >& indices_dists) : radius(radius_), m_indices_dists(indices_dists)
159 inline void clear () { m_indices_dists.clear(); }
161 inline size_t size ()
const {
return m_indices_dists.size(); }
163 inline bool full ()
const {
return true ; }
168 m_indices_dists.push_back(std::make_pair(
index ,dist));
171 inline DistanceType
worstDist ()
const {
return radius; }
186 if (m_indices_dists.empty())
throw std::runtime_error(
"Cannot invoke RadiusResultSet::worst_item() on an empty list of results." );
187 typedef typename std::vector<std::pair<IndexType,DistanceType> >
::const_iterator DistIt;
188 DistIt it = std::max_element(m_indices_dists.begin(), m_indices_dists.end());
197 template <
typename PairType>
198 inline bool operator() (
const PairType &p1,
const PairType &p2)
const {
199 return p1.second < p2.second;
218 fwrite(&
size ,
sizeof (
size_t ), 1, stream);
219 fwrite(&
value [0],
sizeof (T),
size , stream);
226 if (read_cnt !=
count ) {
227 throw std::runtime_error(
"Cannot read from file" );
236 size_t read_cnt = fread(&
size ,
sizeof (
size_t ), 1, stream);
238 throw std::runtime_error(
"Cannot read from file" );
241 read_cnt = fread(&
value [0],
sizeof (T),
size , stream);
242 if (read_cnt!=
size ) {
243 throw std::runtime_error(
"Cannot read from file" );
257 template <
class T,
class DataSource,
typename _DistanceType = T>
265 L1_Adaptor (
const DataSource &_data_source) : data_source(_data_source) { }
270 const T* last =
a +
size ;
271 const T* lastgroup = last - 3;
275 while (
a < lastgroup) {
276 const DistanceType diff0 = std::abs(
a [0] - data_source.kdtree_get_pt(b_idx,d++));
277 const DistanceType diff1 = std::abs(
a [1] - data_source.kdtree_get_pt(b_idx,d++));
278 const DistanceType diff2 = std::abs(
a [2] - data_source.kdtree_get_pt(b_idx,d++));
279 const DistanceType diff3 = std::abs(
a [3] - data_source.kdtree_get_pt(b_idx,d++));
280 result += diff0 + diff1 + diff2 + diff3;
282 if ((worst_dist>0)&&(result>worst_dist)) {
288 result += std::abs( *
a ++ - data_source.kdtree_get_pt(b_idx,d++) );
293 template <
typename U,
typename V>
296 return std::abs(
a -
b );
305 template <
class T,
class DataSource,
typename _DistanceType = T>
313 L2_Adaptor (
const DataSource &_data_source) : data_source(_data_source) { }
318 const T* last =
a +
size ;
319 const T* lastgroup = last - 3;
323 while (
a < lastgroup) {
324 const DistanceType diff0 =
a [0] - data_source.kdtree_get_pt(b_idx,d++);
325 const DistanceType diff1 =
a [1] - data_source.kdtree_get_pt(b_idx,d++);
326 const DistanceType diff2 =
a [2] - data_source.kdtree_get_pt(b_idx,d++);
327 const DistanceType diff3 =
a [3] - data_source.kdtree_get_pt(b_idx,d++);
328 result += diff0 * diff0 + diff1 * diff1 + diff2 * diff2 + diff3 * diff3;
330 if ((worst_dist>0)&&(result>worst_dist)) {
336 const DistanceType diff0 = *
a ++ - data_source.kdtree_get_pt(b_idx,d++);
337 result += diff0 * diff0;
342 template <
typename U,
typename V>
354 template <
class T,
class DataSource,
typename _DistanceType = T>
365 return data_source.kdtree_distance(
a ,b_idx,
size );
368 template <
typename U,
typename V>
377 template <
class T,
class DataSource>
384 template <
class T,
class DataSource>
391 template <
class T,
class DataSource>
406 leaf_max_size(_leaf_max_size)
416 SearchParams (
int checks_IGNORED_ = 32,
float eps_ = 0,
bool sorted_ =
true ) :
417 checks(checks_IGNORED_),
eps (eps_), sorted(sorted_) {}
436 template <
typename T>
439 T* mem =
static_cast< T*
> ( ::malloc(
sizeof (T)*
count ));
503 while (base != NULL) {
504 void *prev = *(
static_cast< void **
> ( base));
526 if (
size > remaining) {
528 wastedMemory += remaining;
535 void * m = ::malloc(blocksize);
537 fprintf (stderr,
"Failed to allocate memory.\n" );
542 static_cast< void **
> (m)[0] = base;
548 remaining = blocksize -
sizeof (
void *) - shift;
549 loc = (
static_cast< char *
> (m) +
sizeof (
void *) + shift);
552 loc =
static_cast< char *
> (loc) +
size ;
567 template <
typename T>
570 T* mem =
static_cast< T*
> (this->malloc(
sizeof (T)*
count ));
606 template <
typename T, std::
size_t N>
628 #if !defined(BOOST_NO_TEMPLATE_PARTIAL_SPECIALIZATION) && !defined(BOOST_MSVC_STD_ITERATOR) && !defined(BOOST_NO_STD_ITERATOR_TRAITS)
631 #elif defined(_MSC_VER) && (_MSC_VER == 1300) && defined(BOOST_DINKUMWARE_STDLIB) && (BOOST_DINKUMWARE_STDLIB == 310)
652 const_reference
at (
size_type i)
const { rangecheck(i);
return elems[i]; }
655 const_reference
front ()
const {
return elems[0]; }
657 const_reference
back ()
const {
return elems[N-1]; }
660 static bool empty () {
return false ; }
664 inline void resize (
const size_t nElements) {
if (nElements!=N)
throw std::logic_error(
"Try to change the size of a CArray." ); }
668 const T*
data ()
const {
return elems; }
679 void assign (
const size_t n ,
const T&
value ) { assert(N==
n );
for (
size_t i=0;i<N;i++) elems[i]=
value ; }
688 template <
int DIM,
typename T>
694 template <
typename T>
739 template <
typename Distance,
class DatasetAdaptor,
int DIM = -1,
typename IndexType =
size_t >
832 dataset(inputData), index_params(
params ), root_node(NULL),
distance (inputData)
834 m_size = dataset.kdtree_get_point_count();
835 m_size_at_index_build =
m_size ;
836 dim = dimensionality;
838 m_leaf_max_size =
params .leaf_max_size;
852 m_size_at_index_build = 0;
862 m_size_at_index_build =
m_size ;
864 computeBoundingBox(root_bbox);
865 root_node = divideTree(0,
m_size , root_bbox );
873 return static_cast< size_t > (DIM>0 ? DIM : dim);
900 template <
typename RESULTSET>
907 throw std::runtime_error(
"[nanoflann] findNeighbors() called before building the index." );
908 float epsError = 1+searchParams.
eps ;
911 dists .assign((DIM>0 ? DIM : dim) ,0);
913 searchLevel(result, vec, root_node, distsq,
dists , epsError);
914 return result.full();
928 resultSet.
init (out_indices, out_distances_sq);
930 return resultSet.
size ();
948 const size_t nFound = radiusSearchCustomCallback(query_point,resultSet,searchParams);
959 template <
class SEARCH_CALLBACK>
962 this->findNeighbors(resultSet, query_point, searchParams);
963 return resultSet.size();
973 m_size = dataset.kdtree_get_point_count();
975 for (
size_t i = 0; i <
m_size ; i++) vind[i] = i;
980 return dataset.kdtree_get_pt(idx,component);
988 save_tree(stream, tree->
child1 );
991 save_tree(stream, tree->
child2 );
1000 if (tree->
child1 !=NULL) {
1001 load_tree(stream, tree->
child1 );
1003 if (tree->
child2 !=NULL) {
1004 load_tree(stream, tree->
child2 );
1011 bbox.
resize ((DIM>0 ? DIM : dim));
1012 if (dataset.kdtree_get_bbox(bbox))
1018 const size_t N = dataset.kdtree_get_point_count();
1019 if (!N)
throw std::runtime_error(
"[nanoflann] computeBoundingBox() called but no data points found." );
1020 for (
int i=0; i<(DIM>0 ? DIM : dim); ++i) {
1022 bbox[i].high = dataset_get(0,i);
1024 for (
size_t k=1; k<N; ++k) {
1025 for (
int i=0; i<(DIM>0 ? DIM : dim); ++i) {
1026 if (dataset_get(k,i)<bbox[i].low) bbox[i].low = dataset_get(k,i);
1027 if (dataset_get(k,i)>bbox[i].high) bbox[i].high = dataset_get(k,i);
1046 if ( (right-left) <=
static_cast< IndexType
> (m_leaf_max_size) ) {
1052 for (
int i=0; i<(DIM>0 ? DIM : dim); ++i) {
1053 bbox[i].low = dataset_get(vind[left],i);
1054 bbox[i].high = dataset_get(vind[left],i);
1056 for (IndexType k=left+1; k<right; ++k) {
1057 for (
int i=0; i<(DIM>0 ? DIM : dim); ++i) {
1058 if (bbox[i].low>dataset_get(vind[k],i)) bbox[i].low=dataset_get(vind[k],i);
1059 if (bbox[i].high<dataset_get(vind[k],i)) bbox[i].high=dataset_get(vind[k],i);
1067 middleSplit_(&vind[0]+left, right-left, idx, cutfeat, cutval, bbox);
1072 left_bbox[cutfeat].high = cutval;
1073 node->
child1 = divideTree(left, left+idx, left_bbox);
1076 right_bbox[cutfeat].low = cutval;
1077 node->
child2 = divideTree(left+idx, right, right_bbox);
1082 for (
int i=0; i<(DIM>0 ? DIM : dim); ++i) {
1083 bbox[i].low =
std::min (left_bbox[i].low, right_bbox[i].low);
1084 bbox[i].high = std::max(left_bbox[i].high, right_bbox[i].high);
1094 min_elem = dataset_get(ind[0],element);
1095 max_elem = dataset_get(ind[0],element);
1096 for (IndexType i=1; i<
count ; ++i) {
1098 if (
val <min_elem) min_elem =
val ;
1099 if (
val >max_elem) max_elem =
val ;
1107 for (
int i=1; i<(DIM>0 ? DIM : dim); ++i) {
1109 if (
span >max_span) {
1115 for (
int i=0; i<(DIM>0 ? DIM : dim); ++i) {
1117 if (
span >(1-EPS)*max_span) {
1119 computeMinMax(ind,
count , cutfeat, min_elem, max_elem);
1121 if (spread>max_spread) {
1123 max_spread = spread;
1128 DistanceType split_val = (bbox[cutfeat].low+bbox[cutfeat].high)/2;
1130 computeMinMax(ind,
count , cutfeat, min_elem, max_elem);
1132 if (split_val<min_elem) cutval = min_elem;
1133 else if (split_val>max_elem) cutval = max_elem;
1134 else cutval = split_val;
1136 IndexType lim1, lim2;
1137 planeSplit(ind,
count , cutfeat, cutval, lim1, lim2);
1158 IndexType right =
count -1;
1160 while (left<=right && dataset_get(ind[left],cutfeat)<cutval) ++left;
1161 while (right && left<=right && dataset_get(ind[right],cutfeat)>=cutval) --right;
1162 if (left>right || !right)
break ;
1163 std::swap(ind[left], ind[right]);
1173 while (left<=right && dataset_get(ind[left],cutfeat)<=cutval) ++left;
1174 while (right && left<=right && dataset_get(ind[right],cutfeat)>cutval) --right;
1175 if (left>right || !right)
break ;
1176 std::swap(ind[left], ind[right]);
1188 for (
int i = 0; i < (DIM>0 ? DIM : dim); ++i) {
1189 if (vec[i] < root_bbox[i].low) {
1190 dists [i] =
distance .accum_dist(vec[i], root_bbox[i].low, i);
1193 if (vec[i] > root_bbox[i].high) {
1194 dists [i] =
distance .accum_dist(vec[i], root_bbox[i].high, i);
1206 template <
class RESULTSET>
1214 for (IndexType i=node->
node_type .
lr .left; i<node->node_type.lr.right; ++i) {
1215 const IndexType
index = vind[i];
1217 if (dist<worst_dist) {
1218 result_set.addPoint(dist,vind[i]);
1233 if ((diff1+diff2)<0) {
1234 bestChild = node->
child1 ;
1235 otherChild = node->
child2 ;
1239 bestChild = node->
child2 ;
1240 otherChild = node->
child1 ;
1245 searchLevel(result_set, vec, bestChild, mindistsq,
dists , epsError);
1248 mindistsq = mindistsq + cut_dist -
dst ;
1249 dists [idx] = cut_dist;
1250 if (mindistsq*epsError<=result_set.worstDist()) {
1251 searchLevel(result_set, vec, otherChild, mindistsq,
dists , epsError);
1268 save_tree(stream, root_node);
1282 load_tree(stream, root_node);
1312 typedef typename Distance::template traits<num_t,self_t>::distance_t
metric_t ;
1321 if (dims!=dimensionality)
throw std::runtime_error(
"Error: 'dimensionality' must match column count in data matrix" );
1322 if (DIM>0 &&
static_cast< int > (dims)!=DIM)
1323 throw std::runtime_error(
"Data set dimensionality does not match the 'DIM' template argument" );
1325 index ->buildIndex();
1343 inline void query (
const num_t *query_point,
const size_t num_closest,
IndexType *out_indices,
num_t *out_distances_sq,
const int = 10)
const
1346 resultSet.
init (out_indices, out_distances_sq);
1362 return m_data_matrix.rows();
1370 const num_t d= p1[i]-m_data_matrix.coeff(idx_p2,i);
1378 return m_data_matrix.coeff(idx,
IndexType (dim));
1384 template <
class BBOX>
int dim
Dimensionality of each data point.
bool findNeighbors(RESULTSET &result, const ElementType *vec, const SearchParams &searchParams) const
An L2-metric KD-tree adaptor for working with data directly stored in an Eigen Matrix,...
Node * child2
Child nodes (both=NULL mean its a leaf node)
const_iterator end() const
const DistanceType radius
num_t kdtree_distance(const num_t *p1, const IndexType idx_p2, IndexType size) const
L2_Simple_Adaptor(const DataSource &_data_source)
DistanceType worstDist() const
reference at(size_type i)
const_reference operator[](size_type i) const
size_t knnSearch(const ElementType *query_point, const size_t num_closest, IndexType *out_indices, DistanceType *out_distances_sq, const int=10) const
Find the "num_closest" nearest neighbors to the query_point[0:dim-1].
EIGEN_STRONG_INLINE iterator begin()
void middleSplit_(IndexType *ind, IndexType count, IndexType &index, int &cutfeat, DistanceType &cutval, const BoundingBox &bbox)
void clear()
Clear the contents of this container.
const Scalar * const_iterator
~KDTreeEigenMatrixAdaptor()
size_t veclen() const
Returns the length of each point in the dataset.
std::pair< IndexType, DistanceType > worst_item() const
Find the worst result (furtherest neighbor) without copying or sorting Pre-conditions: size() > 0.
void addPoint(DistanceType dist, IndexType index)
KNNResultSet(CountType capacity_)
void set_radius_and_clear(const DistanceType r)
Clears the result set and adjusts the search radius.
std::reverse_iterator< iterator > reverse_iterator
operator "<" for std::sort()
KDTreeEigenMatrixAdaptor(const int dimensionality, const MatrixType &mat, const int leaf_max_size=10)
The kd-tree index for the user to call its methods as usual with any other FLANN index.
~PooledAllocator()
Destructor.
PooledAllocator()
Default constructor.
void load_value(FILE *stream, T &value, size_t count=1)
reference operator[](size_type i)
size_t usedMemory() const
Computes the inde memory usage Returns: memory used by the index.
void load_tree(FILE *stream, NodePtr &tree)
void resize(const size_t nElements)
This method has no effects in this class, but raises an exception if the expected size does not match...
void init_vind()
Make sure the auxiliary list vind has the same size than the current dataset, and re-generate if size...
DistanceType operator()(const T *a, const size_t b_idx, size_t size, DistanceType worst_dist=-1) const
const DataSource & data_source
IndexType right
Indices of points in leaf node.
NodePtr divideTree(const IndexType left, const IndexType right, BoundingBox &bbox)
Create a tree node that subdivides the list of vecs from vind[first] to vind[last].
const DataSource & data_source
num_t kdtree_get_pt(const IndexType idx, int dim) const
float eps
search for eps-approximate neighbours (default: 0)
bool kdtree_get_bbox(BBOX &) const
const_iterator begin() const
size_t m_size_at_index_build
Number of points in the dataset when the index was built.
void addPoint(DistanceType dist, IndexType index)
double distance(const TPoint2D &p1, const TPoint2D &p2)
Gets the distance between two points in a 2D space.
T * allocate(size_t count=1)
Allocates (using C's malloc) a generic type T.
void saveIndex(FILE *stream)
Stores the index in a binary file.
Distance::template traits< num_t, self_t >::distance_t metric_t
const_reference front() const
struct nanoflann::KDTreeSingleIndexAdaptor::Node::@23::leaf lr
L1_Adaptor(const DataSource &_data_source)
KDTreeSingleIndexAdaptorParams(size_t _leaf_max_size=10)
ElementType dataset_get(size_t idx, int component) const
Helper accessor to the dataset points:
DistanceType accum_dist(const U a, const V b, int) const
void buildIndex()
Builds the index.
Squared Euclidean (L2) distance functor (suitable for low-dimensionality datasets,...
array_or_vector_selector< DIM, DistanceType >::container_t distance_vector_t
Define "distance_vector_t" as a fixed-size or variable-size container depending on "DIM".
bool sorted
only for radius search, require neighbours sorted by distance (default: true)
_DistanceType DistanceType
RadiusResultSet(DistanceType radius_, std::vector< std::pair< IndexType, DistanceType > > &indices_dists)
GLdouble GLdouble GLdouble r
L2_Simple_Adaptor< T, DataSource > distance_t
Squared Euclidean distance functor (generic version, optimized for high-dimensionality data sets).
size_t size() const
Returns number of points in dataset
T * allocate(const size_t count=1)
Allocates (using this pool) a generic type T.
L1_Adaptor< T, DataSource > distance_t
Search options for KDTreeSingleIndexAdaptor::findNeighbors()
static void rangecheck(size_type i)
const MatrixType & m_data_matrix
const typedef T & const_reference
const_reference at(size_type i) const
bool operator()(const PairType &p1, const PairType &p2) const
PairType will be typically: std::pair<IndexType,DistanceType>
void * malloc(const size_t req_size)
Returns a pointer to a piece of new memory of the given size in bytes allocated from the pool.
const_reference back() const
GLuint GLuint GLsizei count
KDTreeSingleIndexAdaptor(const int dimensionality, const DatasetAdaptor &inputData, const KDTreeSingleIndexAdaptorParams ¶ms=KDTreeSingleIndexAdaptorParams())
KDTree constructor.
Metaprogramming helper traits class for the L2 (Euclidean) metric.
void save_tree(FILE *stream, NodePtr tree)
A result-set class used when performing a radius based search.
const typedef T * const_iterator
std::reverse_iterator< const_iterator > const_reverse_iterator
void init(IndexType *indices_, DistanceType *dists_)
void computeMinMax(IndexType *ind, IndexType count, int element, ElementType &min_elem, ElementType &max_elem)
const DataSource & data_source
size_t radiusSearchCustomCallback(const ElementType *query_point, SEARCH_CALLBACK &resultSet, const SearchParams &searchParams=SearchParams()) const
Just like radiusSearch() but with a custom callback class for each point found in the radius of the q...
L2_Adaptor< T, DataSource > distance_t
std::vector< IndexType > vind
Array of indices to vectors in the dataset.
DistanceType accum_dist(const U a, const V b, int) const
_DistanceType DistanceType
size_t m_size
Number of elements accessed with write access so far.
void query(const num_t *query_point, const size_t num_closest, IndexType *out_indices, num_t *out_distances_sq, const int=10) const
Query for the num_closest closest points to a given point (entered as query_point[0:dim-1]).
_DistanceType DistanceType
struct nanoflann::KDTreeSingleIndexAdaptor::Node::@23::nonleaf sub
DistanceType worstDist() const
void loadIndex(FILE *stream)
Loads a previous index from a binary file.
KDTreeSingleIndexAdaptor< metric_t, self_t, DIM, IndexType > index_t
MatrixType::Index IndexType
std::vector< std::pair< IndexType, DistanceType > > & m_indices_dists
NodePtr root_node
The KD-tree used to find neighbours.
void assign(const T &value)
const self_t & derived() const
Parameters (see README.md)
const KDTreeSingleIndexAdaptorParams index_params
A STL container (as wrapper) for arrays of constant size defined at compile time (class imported from...
GLenum GLenum GLvoid GLvoid GLvoid * span
Metaprogramming helper traits class for the L1 (Manhattan) metric.
const_reverse_iterator rend() const
int fprintf(FILE *fil, const char *format,...) noexcept MRPT_printf_format_check(2
An OS-independent version of fprintf.
union nanoflann::KDTreeSingleIndexAdaptor::Node::@23 node_type
Union used because a node can be either a LEAF node or a non-leaf node, so both data fields are never...
L2_Adaptor(const DataSource &_data_source)
int checks
Ignored parameter (Kept for compatibility with the FLANN interface).
void swap(CArray< T, N > &y)
size_t radiusSearch(const ElementType *query_point, const DistanceType &radius, std::vector< std::pair< IndexType, DistanceType > > &IndicesDists, const SearchParams &searchParams) const
Find all the neighbors to query_point[0:dim-1] within a maximum radius.
GLsizei const GLfloat * value
DistanceType operator()(const T *a, const size_t b_idx, size_t size, DistanceType worst_dist=-1) const
size_t m_size
Number of current poins in the dataset.
~KDTreeSingleIndexAdaptor()
Standard destructor.
const_reverse_iterator rbegin() const
void assign(const size_t n, const T &value)
const DatasetAdaptor & dataset
The dataset used by this index.
KDTreeEigenMatrixAdaptor< MatrixType, DIM, Distance > self_t
std::vector< T > container_t
int divfeat
Dimension used for subdivision.
Metaprogramming helper traits class for the L2_simple (Euclidean) metric.
Manhattan distance functor (generic version, optimized for high-dimensionality data sets).
DistanceType accum_dist(const U a, const V b, int) const
Used to declare fixed-size arrays when DIM>0, dynamically-allocated vectors when DIM=-1.
DistanceType operator()(const T *a, const size_t b_idx, size_t size) const
SearchParams(int checks_IGNORED_=32, float eps_=0, bool sorted_=true)
Note: The first argument (checks_IGNORED_) is ignored, but kept for compatibility with the FLANN inte...
void freeIndex()
Frees the previously-built index.
void searchLevel(RESULTSET &result_set, const ElementType *vec, const NodePtr node, DistanceType mindistsq, distance_vector_t &dists, const float epsError) const
Performs an exact search in the tree starting from a node.
Distance::ElementType ElementType
array_or_vector_selector< DIM, Interval >::container_t BoundingBox
Define "BoundingBox" as a fixed-size or variable-size container depending on "DIM".
CArray< T, DIM > container_t
reverse_iterator rbegin()
void save_value(FILE *stream, const T &value, size_t count=1)
Distance::DistanceType DistanceType
void free_all()
Frees all allocated memory chunks.
PooledAllocator pool
Pooled memory allocator.
DistanceType computeInitialDistances(const ElementType *vec, distance_vector_t &dists) const
GLuint GLuint GLsizei GLenum const GLvoid * indices
void computeBoundingBox(BoundingBox &bbox)
void planeSplit(IndexType *ind, const IndexType count, int cutfeat, DistanceType &cutval, IndexType &lim1, IndexType &lim2)
Subdivide the list of points by a plane perpendicular on axe corresponding to the 'cutfeat' dimension...
const size_t WORDSIZE
Pooled storage allocator.
GLenum const GLfloat * params
GLubyte GLubyte GLubyte a
size_t kdtree_get_point_count() const
static size_type max_size()
std::ptrdiff_t difference_type
Page generated by Doxygen 1.8.17 for MRPT 1.9.9 Git: ad3a9d8ae Tue May 1 23:10:22 2018 -0700 at mié 12 jul 2023 10:03:34 CEST