Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:18:25

0001 /***********************************************************************
0002  * Software License Agreement (BSD License)
0003  *
0004  * Copyright 2008-2009  Marius Muja (mariusm@cs.ubc.ca). All rights reserved.
0005  * Copyright 2008-2009  David G. Lowe (lowe@cs.ubc.ca). All rights reserved.
0006  * Copyright 2011-2016  Jose Luis Blanco (joseluisblancoc@gmail.com).
0007  *   All rights reserved.
0008  *
0009  * THE BSD LICENSE
0010  *
0011  * Redistribution and use in source and binary forms, with or without
0012  * modification, are permitted provided that the following conditions
0013  * are met:
0014  *
0015  * 1. Redistributions of source code must retain the above copyright
0016  *    notice, this list of conditions and the following disclaimer.
0017  * 2. Redistributions in binary form must reproduce the above copyright
0018  *    notice, this list of conditions and the following disclaimer in the
0019  *    documentation and/or other materials provided with the distribution.
0020  *
0021  * THIS SOFTWARE IS PROVIDED BY THE AUTHOR ``AS IS'' AND ANY EXPRESS OR
0022  * IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES
0023  * OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
0024  * IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR ANY DIRECT, INDIRECT,
0025  * INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT
0026  * NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
0027  * DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
0028  * THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
0029  * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF
0030  * THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
0031  *************************************************************************/
0032 
0033 /** \mainpage nanoflann C++ API documentation
0034   *  nanoflann is a C++ header-only library for building KD-Trees, mostly
0035   *  optimized for 2D or 3D point clouds.
0036   *
0037   *  nanoflann does not require compiling or installing, just an
0038   *  #include <nanoflann.hpp> in your code.
0039   *
0040   *  See:
0041   *   - <a href="modules.html" >C++ API organized by modules</a>
0042   *   - <a href="https://github.com/jlblancoc/nanoflann" >Online README</a>
0043   *   - <a href="http://jlblancoc.github.io/nanoflann/" >Doxygen documentation</a>
0044   */
0045 
0046 #ifndef NANOFLANN_HPP_
0047 #define NANOFLANN_HPP_
0048 
0049 #include <algorithm>
0050 #include <cassert>
0051 #include <cmath>    // for abs()
0052 #include <cstdio>   // for fwrite()
0053 #include <cstdlib>  // for abs()
0054 #include <limits>
0055 #include <stdexcept>
0056 #include <vector>
0057 
0058 // Avoid conflicting declaration of min/max macros in windows headers
0059 #if !defined(NOMINMAX) && (defined(_WIN32) || defined(_WIN32_) || defined(WIN32) || defined(_WIN64))
0060 #define NOMINMAX
0061 #ifdef max
0062 #undef max
0063 #undef min
0064 #endif
0065 #endif
0066 
0067 namespace nanoflann
0068 {
0069   /** @addtogroup nanoflann_grp nanoflann C++ library for ANN
0070   *  @{ */
0071 
0072 /** Library version: 0xMmP (M=Major,m=minor,P=patch) */
0073 #define NANOFLANN_VERSION 0x123
0074 
0075   /** @addtogroup result_sets_grp Result set classes
0076       *  @{ */
0077   template <typename DistanceType, typename IndexType = size_t, typename CountType = size_t>
0078   class KNNResultSet
0079   {
0080     IndexType* indices;
0081     DistanceType* dists;
0082     CountType capacity;
0083     CountType count;
0084 
0085    public:
0086     inline KNNResultSet(CountType capacity_)
0087       : indices(0)
0088       , dists(0)
0089       , capacity(capacity_)
0090       , count(0)
0091     {
0092     }
0093 
0094     inline void init(IndexType* indices_, DistanceType* dists_)
0095     {
0096       indices = indices_;
0097       dists = dists_;
0098       count = 0;
0099       if (capacity)
0100         dists[capacity - 1] = (std::numeric_limits<DistanceType>::max)();
0101     }
0102 
0103     inline CountType size() const
0104     {
0105       return count;
0106     }
0107 
0108     inline bool full() const
0109     {
0110       return count == capacity;
0111     }
0112 
0113     inline void addPoint(DistanceType dist, IndexType index)
0114     {
0115       CountType i;
0116       for (i = count; i > 0; --i)
0117       {
0118 #ifdef NANOFLANN_FIRST_MATCH  // If defined and two points have the same distance, the one with the lowest-index will be returned first.
0119         if ((dists[i - 1] > dist) || ((dist == dists[i - 1]) && (indices[i - 1] > index)))
0120         {
0121 #else
0122         if (dists[i - 1] > dist)
0123         {
0124 #endif
0125           if (i < capacity)
0126           {
0127             dists[i] = dists[i - 1];
0128             indices[i] = indices[i - 1];
0129           }
0130         }
0131         else
0132           break;
0133       }
0134       if (i < capacity)
0135       {
0136         dists[i] = dist;
0137         indices[i] = index;
0138       }
0139       if (count < capacity) count++;
0140     }
0141 
0142     inline DistanceType worstDist() const
0143     {
0144       return dists[capacity - 1];
0145     }
0146   };  // namespace nanoflann
0147 
0148   /**
0149      * A result-set class used when performing a radius based search.
0150      */
0151   template <typename DistanceType, typename IndexType = size_t>
0152   class RadiusResultSet
0153   {
0154    public:
0155     const DistanceType radius;
0156 
0157     std::vector<std::pair<IndexType, DistanceType> >& m_indices_dists;
0158 
0159     inline RadiusResultSet(DistanceType radius_, std::vector<std::pair<IndexType, DistanceType> >& indices_dists)
0160       : radius(radius_)
0161       , m_indices_dists(indices_dists)
0162     {
0163       init();
0164     }
0165 
0166     inline ~RadiusResultSet() {}
0167 
0168     inline void init() { clear(); }
0169     inline void clear() { m_indices_dists.clear(); }
0170 
0171     inline size_t size() const { return m_indices_dists.size(); }
0172 
0173     inline bool full() const { return true; }
0174 
0175     inline void addPoint(DistanceType dist, IndexType index)
0176     {
0177       if (dist < radius)
0178         m_indices_dists.push_back(std::make_pair(index, dist));
0179     }
0180 
0181     inline DistanceType worstDist() const { return radius; }
0182 
0183     /** Clears the result set and adjusts the search radius. */
0184     inline void set_radius_and_clear(const DistanceType r)
0185     {
0186       radius = r;
0187       clear();
0188     }
0189 
0190     /**
0191          * Find the worst result (furtherest neighbor) without copying or sorting
0192          * Pre-conditions: size() > 0
0193          */
0194     std::pair<IndexType, DistanceType> worst_item() const
0195     {
0196       if (m_indices_dists.empty()) throw std::runtime_error("Cannot invoke RadiusResultSet::worst_item() on an empty list of results.");
0197       typedef typename std::vector<std::pair<IndexType, DistanceType> >::const_iterator DistIt;
0198       DistIt it = std::max_element(m_indices_dists.begin(), m_indices_dists.end());
0199       return *it;
0200     }
0201   };
0202 
0203   /** operator "<" for std::sort() */
0204   struct IndexDist_Sorter
0205   {
0206     /** PairType will be typically: std::pair<IndexType,DistanceType> */
0207     template <typename PairType>
0208     inline bool operator()(const PairType& p1, const PairType& p2) const
0209     {
0210       return p1.second < p2.second;
0211     }
0212   };
0213 
0214   /** @} */
0215 
0216   /** @addtogroup loadsave_grp Load/save auxiliary functions
0217       * @{ */
0218   template <typename T>
0219   void save_value(FILE* stream, const T& value, size_t count = 1)
0220   {
0221     fwrite(&value, sizeof(value), count, stream);
0222   }
0223 
0224   template <typename T>
0225   void save_value(FILE* stream, const std::vector<T>& value)
0226   {
0227     size_t size = value.size();
0228     fwrite(&size, sizeof(size_t), 1, stream);
0229     fwrite(&value[0], sizeof(T), size, stream);
0230   }
0231 
0232   template <typename T>
0233   void load_value(FILE* stream, T& value, size_t count = 1)
0234   {
0235     size_t read_cnt = fread(&value, sizeof(value), count, stream);
0236     if (read_cnt != count)
0237     {
0238       throw std::runtime_error("Cannot read from file");
0239     }
0240   }
0241 
0242   template <typename T>
0243   void load_value(FILE* stream, std::vector<T>& value)
0244   {
0245     size_t size;
0246     size_t read_cnt = fread(&size, sizeof(size_t), 1, stream);
0247     if (read_cnt != 1)
0248     {
0249       throw std::runtime_error("Cannot read from file");
0250     }
0251     value.resize(size);
0252     read_cnt = fread(&value[0], sizeof(T), size, stream);
0253     if (read_cnt != size)
0254     {
0255       throw std::runtime_error("Cannot read from file");
0256     }
0257   }
0258   /** @} */
0259 
0260   /** @addtogroup metric_grp Metric (distance) classes
0261       * @{ */
0262 
0263   /** Manhattan distance functor (generic version, optimized for high-dimensionality data sets).
0264       *  Corresponding distance traits: nanoflann::metric_L1
0265       * \tparam T Type of the elements (e.g. double, float, uint8_t)
0266       * \tparam _DistanceType Type of distance variables (must be signed) (e.g. float, double, int64_t)
0267       */
0268   template <class T, class DataSource, typename _DistanceType = T>
0269   struct L1_Adaptor
0270   {
0271     typedef T ElementType;
0272     typedef _DistanceType DistanceType;
0273 
0274     const DataSource& data_source;
0275 
0276     L1_Adaptor(const DataSource& _data_source)
0277       : data_source(_data_source)
0278     {
0279     }
0280 
0281     inline DistanceType operator()(const T* a, const size_t b_idx, size_t size, DistanceType worst_dist = -1) const
0282     {
0283       DistanceType result = DistanceType();
0284       const T* last = a + size;
0285       const T* lastgroup = last - 3;
0286       size_t d = 0;
0287 
0288       /* Process 4 items with each loop for efficiency. */
0289       while (a < lastgroup)
0290       {
0291         const DistanceType diff0 = std::abs(a[0] - data_source.kdtree_get_pt(b_idx, d++));
0292         const DistanceType diff1 = std::abs(a[1] - data_source.kdtree_get_pt(b_idx, d++));
0293         const DistanceType diff2 = std::abs(a[2] - data_source.kdtree_get_pt(b_idx, d++));
0294         const DistanceType diff3 = std::abs(a[3] - data_source.kdtree_get_pt(b_idx, d++));
0295         result += diff0 + diff1 + diff2 + diff3;
0296         a += 4;
0297         if ((worst_dist > 0) && (result > worst_dist))
0298         {
0299           return result;
0300         }
0301       }
0302       /* Process last 0-3 components.  Not needed for standard vector lengths. */
0303       while (a < last)
0304       {
0305         result += std::abs(*a++ - data_source.kdtree_get_pt(b_idx, d++));
0306       }
0307       return result;
0308     }
0309 
0310     template <typename U, typename V>
0311     inline DistanceType accum_dist(const U a, const V b, int) const
0312     {
0313       return std::abs(a - b);
0314     }
0315   };
0316 
0317   /** Squared Euclidean distance functor (generic version, optimized for high-dimensionality data sets).
0318       *  Corresponding distance traits: nanoflann::metric_L2
0319       * \tparam T Type of the elements (e.g. double, float, uint8_t)
0320       * \tparam _DistanceType Type of distance variables (must be signed) (e.g. float, double, int64_t)
0321       */
0322   template <class T, class DataSource, typename _DistanceType = T>
0323   struct L2_Adaptor
0324   {
0325     typedef T ElementType;
0326     typedef _DistanceType DistanceType;
0327 
0328     const DataSource& data_source;
0329 
0330     L2_Adaptor(const DataSource& _data_source)
0331       : data_source(_data_source)
0332     {
0333     }
0334 
0335     inline DistanceType operator()(const T* a, const size_t b_idx, size_t size, DistanceType worst_dist = -1) const
0336     {
0337       DistanceType result = DistanceType();
0338       const T* last = a + size;
0339       const T* lastgroup = last - 3;
0340       size_t d = 0;
0341 
0342       /* Process 4 items with each loop for efficiency. */
0343       while (a < lastgroup)
0344       {
0345         const DistanceType diff0 = a[0] - data_source.kdtree_get_pt(b_idx, d++);
0346         const DistanceType diff1 = a[1] - data_source.kdtree_get_pt(b_idx, d++);
0347         const DistanceType diff2 = a[2] - data_source.kdtree_get_pt(b_idx, d++);
0348         const DistanceType diff3 = a[3] - data_source.kdtree_get_pt(b_idx, d++);
0349         result += diff0 * diff0 + diff1 * diff1 + diff2 * diff2 + diff3 * diff3;
0350         a += 4;
0351         if ((worst_dist > 0) && (result > worst_dist))
0352         {
0353           return result;
0354         }
0355       }
0356       /* Process last 0-3 components.  Not needed for standard vector lengths. */
0357       while (a < last)
0358       {
0359         const DistanceType diff0 = *a++ - data_source.kdtree_get_pt(b_idx, d++);
0360         result += diff0 * diff0;
0361       }
0362       return result;
0363     }
0364 
0365     template <typename U, typename V>
0366     inline DistanceType accum_dist(const U a, const V b, int) const
0367     {
0368       return (a - b) * (a - b);
0369     }
0370   };
0371 
0372   /** Squared Euclidean (L2) distance functor (suitable for low-dimensionality datasets, like 2D or 3D point clouds)
0373       *  Corresponding distance traits: nanoflann::metric_L2_Simple
0374       * \tparam T Type of the elements (e.g. double, float, uint8_t)
0375       * \tparam _DistanceType Type of distance variables (must be signed) (e.g. float, double, int64_t)
0376       */
0377   template <class T, class DataSource, typename _DistanceType = T>
0378   struct L2_Simple_Adaptor
0379   {
0380     typedef T ElementType;
0381     typedef _DistanceType DistanceType;
0382 
0383     const DataSource& data_source;
0384 
0385     L2_Simple_Adaptor(const DataSource& _data_source)
0386       : data_source(_data_source)
0387     {
0388     }
0389 
0390     inline DistanceType operator()(const T* a, const size_t b_idx, size_t size) const
0391     {
0392       return data_source.kdtree_distance(a, b_idx, size);
0393     }
0394 
0395     template <typename U, typename V>
0396     inline DistanceType accum_dist(const U a, const V b, int) const
0397     {
0398       return (a - b) * (a - b);
0399     }
0400   };
0401 
0402   /** Metaprogramming helper traits class for the L1 (Manhattan) metric */
0403   struct metric_L1
0404   {
0405     template <class T, class DataSource>
0406     struct traits
0407     {
0408       typedef L1_Adaptor<T, DataSource> distance_t;
0409     };
0410   };
0411   /** Metaprogramming helper traits class for the L2 (Euclidean) metric */
0412   struct metric_L2
0413   {
0414     template <class T, class DataSource>
0415     struct traits
0416     {
0417       typedef L2_Adaptor<T, DataSource> distance_t;
0418     };
0419   };
0420   /** Metaprogramming helper traits class for the L2_simple (Euclidean) metric */
0421   struct metric_L2_Simple
0422   {
0423     template <class T, class DataSource>
0424     struct traits
0425     {
0426       typedef L2_Simple_Adaptor<T, DataSource> distance_t;
0427     };
0428   };
0429 
0430   /** @} */
0431 
0432   /** @addtogroup param_grp Parameter structs
0433       * @{ */
0434 
0435   /**  Parameters (see README.md) */
0436   struct KDTreeSingleIndexAdaptorParams
0437   {
0438     KDTreeSingleIndexAdaptorParams(size_t _leaf_max_size = 10)
0439       : leaf_max_size(_leaf_max_size)
0440     {
0441     }
0442 
0443     size_t leaf_max_size;
0444   };
0445 
0446   /** Search options for KDTreeSingleIndexAdaptor::findNeighbors() */
0447   struct SearchParams
0448   {
0449     /** Note: The first argument (checks_IGNORED_) is ignored, but kept for compatibility with the FLANN interface */
0450     SearchParams(int checks_IGNORED_ = 32, float eps_ = 0, bool sorted_ = true)
0451       : checks(checks_IGNORED_)
0452       , eps(eps_)
0453       , sorted(sorted_)
0454     {
0455     }
0456 
0457     int checks;   //!< Ignored parameter (Kept for compatibility with the FLANN interface).
0458     float eps;    //!< search for eps-approximate neighbours (default: 0)
0459     bool sorted;  //!< only for radius search, require neighbours sorted by distance (default: true)
0460   };
0461   /** @} */
0462 
0463   /** @addtogroup memalloc_grp Memory allocation
0464       * @{ */
0465 
0466   /**
0467      * Allocates (using C's malloc) a generic type T.
0468      *
0469      * Params:
0470      *     count = number of instances to allocate.
0471      * Returns: pointer (of type T*) to memory buffer
0472      */
0473   template <typename T>
0474   inline T* allocate(size_t count = 1)
0475   {
0476     T* mem = static_cast<T*>(::malloc(sizeof(T) * count));
0477     return mem;
0478   }
0479 
0480   /**
0481      * Pooled storage allocator
0482      *
0483      * The following routines allow for the efficient allocation of storage in
0484      * small chunks from a specified pool.  Rather than allowing each structure
0485      * to be freed individually, an entire pool of storage is freed at once.
0486      * This method has two advantages over just using malloc() and free().  First,
0487      * it is far more efficient for allocating small objects, as there is
0488      * no overhead for remembering all the information needed to free each
0489      * object or consolidating fragmented memory.  Second, the decision about
0490      * how long to keep an object is made at the time of allocation, and there
0491      * is no need to track down all the objects to free them.
0492      *
0493      */
0494 
0495   const size_t WORDSIZE = 16;
0496   const size_t BLOCKSIZE = 8192;
0497 
0498   class PooledAllocator
0499   {
0500     /* We maintain memory alignment to word boundaries by requiring that all
0501             allocations be in multiples of the machine wordsize.  */
0502     /* Size of machine word in bytes.  Must be power of 2. */
0503     /* Minimum number of bytes requested at a time from the system.  Must be multiple of WORDSIZE. */
0504 
0505     size_t remaining; /* Number of bytes left in current block of storage. */
0506     void* base;       /* Pointer to base of current block of storage. */
0507     void* loc;        /* Current location in block to next allocate memory. */
0508 
0509     void internal_init()
0510     {
0511       remaining = 0;
0512       base = NULL;
0513       usedMemory = 0;
0514       wastedMemory = 0;
0515     }
0516 
0517    public:
0518     size_t usedMemory;
0519     size_t wastedMemory;
0520 
0521     /**
0522             Default constructor. Initializes a new pool.
0523          */
0524 // cppcheck-suppress uninitMemberVar
0525     PooledAllocator()
0526     {
0527       internal_init();
0528     }
0529 
0530     /**
0531          * Destructor. Frees all the memory allocated in this pool.
0532          */
0533     ~PooledAllocator()
0534     {
0535       free_all();
0536     }
0537 
0538     /** Frees all allocated memory chunks */
0539     void free_all()
0540     {
0541       while (base != NULL)
0542       {
0543         void* prev = *(static_cast<void**>(base)); /* Get pointer to prev block. */
0544         ::free(base);
0545         base = prev;
0546       }
0547       internal_init();
0548     }
0549 
0550     /**
0551          * Returns a pointer to a piece of new memory of the given size in bytes
0552          * allocated from the pool.
0553          */
0554     void* malloc(const size_t req_size)
0555     {
0556       /* Round size up to a multiple of wordsize.  The following expression
0557                 only works for WORDSIZE that is a power of 2, by masking last bits of
0558                 incremented size to zero.
0559              */
0560       const size_t size = (req_size + (WORDSIZE - 1)) & ~(WORDSIZE - 1);
0561 
0562       /* Check whether a new block must be allocated.  Note that the first word
0563                 of a block is reserved for a pointer to the previous block.
0564              */
0565       if (size > remaining)
0566       {
0567         wastedMemory += remaining;
0568 
0569         /* Allocate new storage. */
0570         const size_t blocksize = (size + sizeof(void*) + (WORDSIZE - 1) > BLOCKSIZE) ? size + sizeof(void*) + (WORDSIZE - 1) : BLOCKSIZE;
0571 
0572         // use the standard C malloc to allocate memory
0573         void* m = ::malloc(blocksize);
0574         if (!m)
0575         {
0576           fprintf(stderr, "Failed to allocate memory.\n");
0577           return NULL;
0578         }
0579 
0580         /* Fill first word of new block with pointer to previous block. */
0581         static_cast<void**>(m)[0] = base;
0582         base = m;
0583 
0584         size_t shift = 0;
0585         //int size_t = (WORDSIZE - ( (((size_t)m) + sizeof(void*)) & (WORDSIZE-1))) & (WORDSIZE-1);
0586 
0587         remaining = blocksize - sizeof(void*) - shift;
0588         loc = (static_cast<char*>(m) + sizeof(void*) + shift);
0589       }
0590       void* rloc = loc;
0591       loc = static_cast<char*>(loc) + size;
0592       remaining -= size;
0593 
0594       usedMemory += size;
0595 
0596       return rloc;
0597     }
0598 
0599     /**
0600          * Allocates (using this pool) a generic type T.
0601          *
0602          * Params:
0603          *     count = number of instances to allocate.
0604          * Returns: pointer (of type T*) to memory buffer
0605          */
0606     template <typename T>
0607     T* allocate(const size_t count = 1)
0608     {
0609       T* mem = static_cast<T*>(this->malloc(sizeof(T) * count));
0610       return mem;
0611     }
0612   };
0613   /** @} */
0614 
0615   /** @addtogroup nanoflann_metaprog_grp Auxiliary metaprogramming stuff
0616       * @{ */
0617 
0618   // ----------------  CArray -------------------------
0619   /** A STL container (as wrapper) for arrays of constant size defined at compile time (class imported from the MRPT project)
0620      * This code is an adapted version from Boost, modifed for its integration
0621      *  within MRPT (JLBC, Dec/2009) (Renamed array -> CArray to avoid possible potential conflicts).
0622      * See
0623      *      http://www.josuttis.com/cppcode
0624      * for details and the latest version.
0625      * See
0626      *      http://www.boost.org/libs/array for Documentation.
0627      * for documentation.
0628      *
0629      * (C) Copyright Nicolai M. Josuttis 2001.
0630      * Permission to copy, use, modify, sell and distribute this software
0631      * is granted provided this copyright notice appears in all copies.
0632      * This software is provided "as is" without express or implied
0633      * warranty, and with no claim as to its suitability for any purpose.
0634      *
0635      * 29 Jan 2004 - minor fixes (Nico Josuttis)
0636      * 04 Dec 2003 - update to synch with library TR1 (Alisdair Meredith)
0637      * 23 Aug 2002 - fix for Non-MSVC compilers combined with MSVC libraries.
0638      * 05 Aug 2001 - minor update (Nico Josuttis)
0639      * 20 Jan 2001 - STLport fix (Beman Dawes)
0640      * 29 Sep 2000 - Initial Revision (Nico Josuttis)
0641      *
0642      * Jan 30, 2004
0643      */
0644   template <typename T, std::size_t N>
0645   class CArray
0646   {
0647    public:
0648     T elems[N];  // fixed-size array of elements of type T
0649 
0650    public:
0651     // type definitions
0652     typedef T value_type;
0653     typedef T* iterator;
0654     typedef const T* const_iterator;
0655     typedef T& reference;
0656     typedef const T& const_reference;
0657     typedef std::size_t size_type;
0658     typedef std::ptrdiff_t difference_type;
0659 
0660     // iterator support
0661     inline iterator begin() { return elems; }
0662     inline const_iterator begin() const { return elems; }
0663     inline iterator end() { return elems + N; }
0664     inline const_iterator end() const { return elems + N; }
0665 
0666     // reverse iterator support
0667 #if !defined(BOOST_NO_TEMPLATE_PARTIAL_SPECIALIZATION) && !defined(BOOST_MSVC_STD_ITERATOR) && !defined(BOOST_NO_STD_ITERATOR_TRAITS)
0668     typedef std::reverse_iterator<iterator> reverse_iterator;
0669     typedef std::reverse_iterator<const_iterator> const_reverse_iterator;
0670 #elif defined(_MSC_VER) && (_MSC_VER == 1300) && defined(BOOST_DINKUMWARE_STDLIB) && (BOOST_DINKUMWARE_STDLIB == 310)
0671     // workaround for broken reverse_iterator in VC7
0672     typedef std::reverse_iterator<std::_Ptrit<value_type, difference_type, iterator,
0673                                               reference, iterator, reference> >
0674         reverse_iterator;
0675     typedef std::reverse_iterator<std::_Ptrit<value_type, difference_type, const_iterator,
0676                                               const_reference, iterator, reference> >
0677         const_reverse_iterator;
0678 #else
0679   // workaround for broken reverse_iterator implementations
0680   typedef std::reverse_iterator<iterator, T> reverse_iterator;
0681   typedef std::reverse_iterator<const_iterator, T> const_reverse_iterator;
0682 #endif
0683 
0684     reverse_iterator rbegin()
0685     {
0686       return reverse_iterator(end());
0687     }
0688     const_reverse_iterator rbegin() const { return const_reverse_iterator(end()); }
0689     reverse_iterator rend() { return reverse_iterator(begin()); }
0690     const_reverse_iterator rend() const { return const_reverse_iterator(begin()); }
0691     // operator[]
0692     inline reference operator[](size_type i) { return elems[i]; }
0693     inline const_reference operator[](size_type i) const { return elems[i]; }
0694     // at() with range check
0695     reference at(size_type i)
0696     {
0697       rangecheck(i);
0698       return elems[i];
0699     }
0700     const_reference at(size_type i) const
0701     {
0702       rangecheck(i);
0703       return elems[i];
0704     }
0705     // front() and back()
0706     reference front() { return elems[0]; }
0707     const_reference front() const { return elems[0]; }
0708     reference back() { return elems[N - 1]; }
0709     const_reference back() const { return elems[N - 1]; }
0710     // size is constant
0711     static inline size_type size() { return N; }
0712     static bool empty() { return false; }
0713     static size_type max_size() { return N; }
0714     enum
0715     {
0716       static_size = N
0717     };
0718     /** This method has no effects in this class, but raises an exception if the expected size does not match */
0719     inline void resize(const size_t nElements)
0720     {
0721       if (nElements != N) throw std::logic_error("Try to change the size of a CArray.");
0722     }
0723     // swap (note: linear complexity in N, constant for given instantiation)
0724     void swap(CArray<T, N>& y) { std::swap_ranges(begin(), end(), y.begin()); }
0725     // direct access to data (read-only)
0726     const T* data() const { return elems; }
0727     // use array as C array (direct read/write access to data)
0728     T* data() { return elems; }
0729     // assignment with type conversion
0730     template <typename T2>
0731 // cppcheck-suppress *
0732     CArray<T, N>& operator=(const CArray<T2, N>& rhs)
0733     {
0734       std::copy(rhs.begin(), rhs.end(), begin());
0735       return *this;
0736     }
0737     // assign one value to all elements
0738     inline void assign(const T& value)
0739     {
0740       for (size_t i = 0; i < N; i++) elems[i] = value;
0741     }
0742     // assign (compatible with std::vector's one) (by JLBC for MRPT)
0743     void assign(const size_t n, const T& value)
0744     {
0745       assert(N == n);
0746       for (size_t i = 0; i < N; i++) elems[i] = value;
0747     }
0748 
0749    private:
0750     // check range (may be private because it is static)
0751     static void rangecheck(size_type i)
0752     {
0753       if (i >= size())
0754       {
0755         throw std::out_of_range("CArray<>: index out of range");
0756       }
0757     }
0758   };  // end of CArray
0759 
0760   /** Used to declare fixed-size arrays when DIM>0, dynamically-allocated vectors when DIM=-1.
0761       * Fixed size version for a generic DIM:
0762       */
0763   template <int DIM, typename T>
0764   struct array_or_vector_selector
0765   {
0766     typedef CArray<T, DIM> container_t;
0767   };
0768   /** Dynamic size version */
0769   template <typename T>
0770   struct array_or_vector_selector<-1, T>
0771   {
0772     typedef std::vector<T> container_t;
0773   };
0774   /** @} */
0775 
0776   /** @addtogroup kdtrees_grp KD-tree classes and adaptors
0777       * @{ */
0778 
0779   /** kd-tree index
0780      *
0781      * Contains the k-d trees and other information for indexing a set of points
0782      * for nearest-neighbor matching.
0783      *
0784      *  The class "DatasetAdaptor" must provide the following interface (can be non-virtual, inlined methods):
0785      *
0786      *  \code
0787      *   // Must return the number of data poins
0788      *   inline size_t kdtree_get_point_count() const { ... }
0789      *
0790      *   // [Only if using the metric_L2_Simple type] Must return the Euclidean (L2) distance between the vector "p1[0:size-1]" and the data point with index "idx_p2" stored in the class:
0791      *   inline DistanceType kdtree_distance(const T *p1, const size_t idx_p2,size_t size) const { ... }
0792      *
0793      *   // Must return the dim'th component of the idx'th point in the class:
0794      *   inline T kdtree_get_pt(const size_t idx, int dim) const { ... }
0795      *
0796      *   // Optional bounding-box computation: return false to default to a standard bbox computation loop.
0797      *   //   Return true if the BBOX was already computed by the class and returned in "bb" so it can be avoided to redo it again.
0798      *   //   Look at bb.size() to find out the expected dimensionality (e.g. 2 or 3 for point clouds)
0799      *   template <class BBOX>
0800      *   bool kdtree_get_bbox(BBOX &bb) const
0801      *   {
0802      *      bb[0].low = ...; bb[0].high = ...;  // 0th dimension limits
0803      *      bb[1].low = ...; bb[1].high = ...;  // 1st dimension limits
0804      *      ...
0805      *      return true;
0806      *   }
0807      *
0808      *  \endcode
0809      *
0810      * \tparam DatasetAdaptor The user-provided adaptor (see comments above).
0811      * \tparam Distance The distance metric to use: nanoflann::metric_L1, nanoflann::metric_L2, nanoflann::metric_L2_Simple, etc.
0812      * \tparam DIM Dimensionality of data points (e.g. 3 for 3D points)
0813      * \tparam IndexType Will be typically size_t or int
0814      */
0815   template <typename Distance, class DatasetAdaptor, int DIM = -1, typename IndexType = size_t>
0816   class KDTreeSingleIndexAdaptor
0817   {
0818    private:
0819     /** Hidden copy constructor, to disallow copying indices (Not implemented) */
0820     KDTreeSingleIndexAdaptor(const KDTreeSingleIndexAdaptor<Distance, DatasetAdaptor, DIM, IndexType>&);
0821 
0822    public:
0823     typedef typename Distance::ElementType ElementType;
0824     typedef typename Distance::DistanceType DistanceType;
0825 
0826    protected:
0827     /**
0828          *  Array of indices to vectors in the dataset.
0829          */
0830     std::vector<IndexType> vind;
0831 
0832     size_t m_leaf_max_size;
0833 
0834     /**
0835          * The dataset used by this index
0836          */
0837     const DatasetAdaptor& dataset;  //!< The source of our data
0838 
0839     const KDTreeSingleIndexAdaptorParams index_params;
0840 
0841     size_t m_size;                 //!< Number of current poins in the dataset
0842     size_t m_size_at_index_build;  //!< Number of points in the dataset when the index was built
0843     int dim;                       //!< Dimensionality of each data point
0844 
0845     /*--------------------- Internal Data Structures --------------------------*/
0846     struct Node
0847     {
0848       /** Union used because a node can be either a LEAF node or a non-leaf node, so both data fields are never used simultaneously */
0849       union {
0850         struct leaf
0851         {
0852           IndexType left, right;  //!< Indices of points in leaf node
0853         } lr;
0854         struct nonleaf
0855         {
0856           int divfeat;                   //!< Dimension used for subdivision.
0857           DistanceType divlow, divhigh;  //!< The values used for subdivision.
0858         } sub;
0859       } node_type;
0860       Node *child1, *child2;  //!< Child nodes (both=NULL mean its a leaf node)
0861     };
0862     typedef Node* NodePtr;
0863 
0864     struct Interval
0865     {
0866       ElementType low, high;
0867     };
0868 
0869     /** Define "BoundingBox" as a fixed-size or variable-size container depending on "DIM" */
0870     typedef typename array_or_vector_selector<DIM, Interval>::container_t BoundingBox;
0871 
0872     /** Define "distance_vector_t" as a fixed-size or variable-size container depending on "DIM" */
0873     typedef typename array_or_vector_selector<DIM, DistanceType>::container_t distance_vector_t;
0874 
0875     /** The KD-tree used to find neighbours */
0876     NodePtr root_node;
0877     BoundingBox root_bbox;
0878 
0879     /**
0880          * Pooled memory allocator.
0881          *
0882          * Using a pooled memory allocator is more efficient
0883          * than allocating memory directly when there is a large
0884          * number small of memory allocations.
0885          */
0886     PooledAllocator pool;
0887 
0888    public:
0889     Distance distance;
0890 
0891     /**
0892          * KDTree constructor
0893          *
0894          * Refer to docs in README.md or online in https://github.com/jlblancoc/nanoflann
0895          *
0896          * The KD-Tree point dimension (the length of each point in the datase, e.g. 3 for 3D points)
0897          * is determined by means of:
0898          *  - The \a DIM template parameter if >0 (highest priority)
0899          *  - Otherwise, the \a dimensionality parameter of this constructor.
0900          *
0901          * @param inputData Dataset with the input features
0902          * @param params Basically, the maximum leaf node size
0903          */
0904     KDTreeSingleIndexAdaptor(const int dimensionality, const DatasetAdaptor& inputData, const KDTreeSingleIndexAdaptorParams& params = KDTreeSingleIndexAdaptorParams())
0905       : dataset(inputData)
0906       , index_params(params)
0907       , root_node(NULL)
0908       , distance(inputData)
0909     {
0910       m_size = dataset.kdtree_get_point_count();
0911       m_size_at_index_build = m_size;
0912       dim = dimensionality;
0913       if (DIM > 0) dim = DIM;
0914       m_leaf_max_size = params.leaf_max_size;
0915 
0916       // Create a permutable array of indices to the input vectors.
0917       init_vind();
0918     }
0919 
0920     /** Standard destructor */
0921     ~KDTreeSingleIndexAdaptor() {}
0922 
0923     /** Frees the previously-built index. Automatically called within buildIndex(). */
0924     void freeIndex()
0925     {
0926       pool.free_all();
0927       root_node = NULL;
0928       m_size_at_index_build = 0;
0929     }
0930 
0931     /**
0932          * Builds the index
0933          */
0934     void buildIndex()
0935     {
0936       init_vind();
0937       freeIndex();
0938       m_size_at_index_build = m_size;
0939       if (m_size == 0) return;
0940       computeBoundingBox(root_bbox);
0941       root_node = divideTree(0, m_size, root_bbox);  // construct the tree
0942     }
0943 
0944     /** Returns number of points in dataset  */
0945     size_t size() const { return m_size; }
0946 
0947     /** Returns the length of each point in the dataset */
0948     size_t veclen() const
0949     {
0950       return static_cast<size_t>(DIM > 0 ? DIM : dim);
0951     }
0952 
0953     /**
0954          * Computes the inde memory usage
0955          * Returns: memory used by the index
0956          */
0957     size_t usedMemory() const
0958     {
0959       return pool.usedMemory + pool.wastedMemory + dataset.kdtree_get_point_count() * sizeof(IndexType);  // pool memory and vind array memory
0960     }
0961 
0962     /** \name Query methods
0963           * @{ */
0964 
0965     /**
0966          * Find set of nearest neighbors to vec[0:dim-1]. Their indices are stored inside
0967          * the result object.
0968          *
0969          * Params:
0970          *     result = the result object in which the indices of the nearest-neighbors are stored
0971          *     vec = the vector for which to search the nearest neighbors
0972          *
0973          * \tparam RESULTSET Should be any ResultSet<DistanceType>
0974          * \return  True if the requested neighbors could be found.
0975          * \sa knnSearch, radiusSearch
0976          */
0977     template <typename RESULTSET>
0978     bool findNeighbors(RESULTSET& result, const ElementType* vec, const SearchParams& searchParams) const
0979     {
0980       assert(vec);
0981       if (size() == 0)
0982         return false;
0983       if (!root_node)
0984         throw std::runtime_error("[nanoflann] findNeighbors() called before building the index.");
0985       float epsError = 1 + searchParams.eps;
0986 
0987       distance_vector_t dists;                 // fixed or variable-sized container (depending on DIM)
0988       dists.assign((DIM > 0 ? DIM : dim), 0);  // Fill it with zeros.
0989       DistanceType distsq = computeInitialDistances(vec, dists);
0990       searchLevel(result, vec, root_node, distsq, dists, epsError);  // "count_leaf" parameter removed since was neither used nor returned to the user.
0991       return result.full();
0992     }
0993 
0994     /**
0995          * Find the "num_closest" nearest neighbors to the \a query_point[0:dim-1]. Their indices are stored inside
0996          * the result object.
0997          *  \sa radiusSearch, findNeighbors
0998          * \note nChecks_IGNORED is ignored but kept for compatibility with the original FLANN interface.
0999          * \return Number `N` of valid points in the result set. Only the first `N` entries in `out_indices` and `out_distances_sq` will be valid. 
1000          *         Return may be less than `num_closest` only if the number of elements in the tree is less than `num_closest`.
1001          */
1002     size_t knnSearch(const ElementType* query_point, const size_t num_closest, IndexType* out_indices, DistanceType* out_distances_sq, const int /* nChecks_IGNORED */ = 10) const
1003     {
1004       nanoflann::KNNResultSet<DistanceType, IndexType> resultSet(num_closest);
1005       resultSet.init(out_indices, out_distances_sq);
1006       this->findNeighbors(resultSet, query_point, nanoflann::SearchParams());
1007       return resultSet.size();
1008     }
1009 
1010     /**
1011          * Find all the neighbors to \a query_point[0:dim-1] within a maximum radius.
1012          *  The output is given as a vector of pairs, of which the first element is a point index and the second the corresponding distance.
1013          *  Previous contents of \a IndicesDists are cleared.
1014          *
1015          *  If searchParams.sorted==true, the output list is sorted by ascending distances.
1016          *
1017          *  For a better performance, it is advisable to do a .reserve() on the vector if you have any wild guess about the number of expected matches.
1018          *
1019          *  \sa knnSearch, findNeighbors, radiusSearchCustomCallback
1020          * \return The number of points within the given radius (i.e. indices.size() or dists.size() )
1021          */
1022     size_t radiusSearch(const ElementType* query_point, const DistanceType& radius, std::vector<std::pair<IndexType, DistanceType> >& IndicesDists, const SearchParams& searchParams) const
1023     {
1024       RadiusResultSet<DistanceType, IndexType> resultSet(radius, IndicesDists);
1025       const size_t nFound = radiusSearchCustomCallback(query_point, resultSet, searchParams);
1026       if (searchParams.sorted)
1027         std::sort(IndicesDists.begin(), IndicesDists.end(), IndexDist_Sorter());
1028       return nFound;
1029     }
1030 
1031     /**
1032          * Just like radiusSearch() but with a custom callback class for each point found in the radius of the query.
1033          * See the source of RadiusResultSet<> as a start point for your own classes.
1034          * \sa radiusSearch
1035          */
1036     template <class SEARCH_CALLBACK>
1037     size_t radiusSearchCustomCallback(const ElementType* query_point, SEARCH_CALLBACK& resultSet, const SearchParams& searchParams = SearchParams()) const
1038     {
1039       this->findNeighbors(resultSet, query_point, searchParams);
1040       return resultSet.size();
1041     }
1042 
1043     /** @} */
1044 
1045    private:
1046     /** Make sure the auxiliary list \a vind has the same size than the current dataset, and re-generate if size has changed. */
1047     void init_vind()
1048     {
1049       // Create a permutable array of indices to the input vectors.
1050       m_size = dataset.kdtree_get_point_count();
1051       if (vind.size() != m_size) vind.resize(m_size);
1052       for (size_t i = 0; i < m_size; i++) vind[i] = i;
1053     }
1054 
1055     /// Helper accessor to the dataset points:
1056     inline ElementType dataset_get(size_t idx, int component) const
1057     {
1058       return dataset.kdtree_get_pt(idx, component);
1059     }
1060 
1061     void save_tree(FILE* stream, NodePtr tree)
1062     {
1063       save_value(stream, *tree);
1064       if (tree->child1 != NULL)
1065       {
1066         save_tree(stream, tree->child1);
1067       }
1068       if (tree->child2 != NULL)
1069       {
1070         save_tree(stream, tree->child2);
1071       }
1072     }
1073 
1074     void load_tree(FILE* stream, NodePtr& tree)
1075     {
1076       tree = pool.allocate<Node>();
1077       load_value(stream, *tree);
1078       if (tree->child1 != NULL)
1079       {
1080         load_tree(stream, tree->child1);
1081       }
1082       if (tree->child2 != NULL)
1083       {
1084         load_tree(stream, tree->child2);
1085       }
1086     }
1087 
1088     void computeBoundingBox(BoundingBox& bbox)
1089     {
1090       bbox.resize((DIM > 0 ? DIM : dim));
1091       if (dataset.kdtree_get_bbox(bbox))
1092       {
1093         // Done! It was implemented in derived class
1094       }
1095       else
1096       {
1097         const size_t N = dataset.kdtree_get_point_count();
1098         if (!N) throw std::runtime_error("[nanoflann] computeBoundingBox() called but no data points found.");
1099         for (int i = 0; i < (DIM > 0 ? DIM : dim); ++i)
1100         {
1101           bbox[i].low =
1102               bbox[i].high = dataset_get(0, i);
1103         }
1104         for (size_t k = 1; k < N; ++k)
1105         {
1106           for (int i = 0; i < (DIM > 0 ? DIM : dim); ++i)
1107           {
1108             if (dataset_get(k, i) < bbox[i].low) bbox[i].low = dataset_get(k, i);
1109             if (dataset_get(k, i) > bbox[i].high) bbox[i].high = dataset_get(k, i);
1110           }
1111         }
1112       }
1113     }
1114 
1115     /**
1116          * Create a tree node that subdivides the list of vecs from vind[first]
1117          * to vind[last].  The routine is called recursively on each sublist.
1118          *
1119          * @param left index of the first vector
1120          * @param right index of the last vector
1121          */
1122     NodePtr divideTree(const IndexType left, const IndexType right, BoundingBox& bbox)
1123     {
1124       NodePtr node = pool.allocate<Node>();  // allocate memory
1125 
1126       /* If too few exemplars remain, then make this a leaf node. */
1127       if ((right - left) <= static_cast<IndexType>(m_leaf_max_size))
1128       {
1129         node->child1 = node->child2 = NULL; /* Mark as leaf node. */
1130         node->node_type.lr.left = left;
1131         node->node_type.lr.right = right;
1132 
1133         // compute bounding-box of leaf points
1134         for (int i = 0; i < (DIM > 0 ? DIM : dim); ++i)
1135         {
1136           bbox[i].low = dataset_get(vind[left], i);
1137           bbox[i].high = dataset_get(vind[left], i);
1138         }
1139         for (IndexType k = left + 1; k < right; ++k)
1140         {
1141           for (int i = 0; i < (DIM > 0 ? DIM : dim); ++i)
1142           {
1143             if (bbox[i].low > dataset_get(vind[k], i)) bbox[i].low = dataset_get(vind[k], i);
1144             if (bbox[i].high < dataset_get(vind[k], i)) bbox[i].high = dataset_get(vind[k], i);
1145           }
1146         }
1147       }
1148       else
1149       {
1150         IndexType idx;
1151         int cutfeat;
1152         DistanceType cutval;
1153         middleSplit_(&vind[0] + left, right - left, idx, cutfeat, cutval, bbox);
1154 
1155         node->node_type.sub.divfeat = cutfeat;
1156 
1157         BoundingBox left_bbox(bbox);
1158         left_bbox[cutfeat].high = cutval;
1159         node->child1 = divideTree(left, left + idx, left_bbox);
1160 
1161         BoundingBox right_bbox(bbox);
1162         right_bbox[cutfeat].low = cutval;
1163         node->child2 = divideTree(left + idx, right, right_bbox);
1164 
1165         node->node_type.sub.divlow = left_bbox[cutfeat].high;
1166         node->node_type.sub.divhigh = right_bbox[cutfeat].low;
1167 
1168         for (int i = 0; i < (DIM > 0 ? DIM : dim); ++i)
1169         {
1170           bbox[i].low = std::min(left_bbox[i].low, right_bbox[i].low);
1171           bbox[i].high = std::max(left_bbox[i].high, right_bbox[i].high);
1172         }
1173       }
1174 
1175       return node;
1176     }
1177 
1178     void computeMinMax(IndexType* ind, IndexType count, int element, ElementType& min_elem, ElementType& max_elem)
1179     {
1180       min_elem = dataset_get(ind[0], element);
1181       max_elem = dataset_get(ind[0], element);
1182       for (IndexType i = 1; i < count; ++i)
1183       {
1184         ElementType val = dataset_get(ind[i], element);
1185         if (val < min_elem) min_elem = val;
1186         if (val > max_elem) max_elem = val;
1187       }
1188     }
1189 
1190     void middleSplit_(IndexType* ind, IndexType count, IndexType& index, int& cutfeat, DistanceType& cutval, const BoundingBox& bbox)
1191     {
1192       const DistanceType EPS = static_cast<DistanceType>(0.00001);
1193       ElementType max_span = bbox[0].high - bbox[0].low;
1194       for (int i = 1; i < (DIM > 0 ? DIM : dim); ++i)
1195       {
1196         ElementType span = bbox[i].high - bbox[i].low;
1197         if (span > max_span)
1198         {
1199           max_span = span;
1200         }
1201       }
1202       ElementType max_spread = -1;
1203       cutfeat = 0;
1204       for (int i = 0; i < (DIM > 0 ? DIM : dim); ++i)
1205       {
1206         ElementType span = bbox[i].high - bbox[i].low;
1207         if (span > (1 - EPS) * max_span)
1208         {
1209           ElementType min_elem, max_elem;
1210           computeMinMax(ind, count, i, min_elem, max_elem);
1211           ElementType spread = max_elem - min_elem;
1212           ;
1213           if (spread > max_spread)
1214           {
1215             cutfeat = i;
1216             max_spread = spread;
1217           }
1218         }
1219       }
1220       // split in the middle
1221       DistanceType split_val = (bbox[cutfeat].low + bbox[cutfeat].high) / 2;
1222       ElementType min_elem, max_elem;
1223       computeMinMax(ind, count, cutfeat, min_elem, max_elem);
1224 
1225       if (split_val < min_elem)
1226         cutval = min_elem;
1227       else if (split_val > max_elem)
1228         cutval = max_elem;
1229       else
1230         cutval = split_val;
1231 
1232       IndexType lim1, lim2;
1233       planeSplit(ind, count, cutfeat, cutval, lim1, lim2);
1234 
1235       if (lim1 > count / 2)
1236         index = lim1;
1237       else if (lim2 < count / 2)
1238         index = lim2;
1239       else
1240         index = count / 2;
1241     }
1242 
1243     /**
1244          *  Subdivide the list of points by a plane perpendicular on axe corresponding
1245          *  to the 'cutfeat' dimension at 'cutval' position.
1246          *
1247          *  On return:
1248          *  dataset[ind[0..lim1-1]][cutfeat]<cutval
1249          *  dataset[ind[lim1..lim2-1]][cutfeat]==cutval
1250          *  dataset[ind[lim2..count]][cutfeat]>cutval
1251          */
1252     void planeSplit(IndexType* ind, const IndexType count, int cutfeat, DistanceType& cutval, IndexType& lim1, IndexType& lim2)
1253     {
1254       /* Move vector indices for left subtree to front of list. */
1255       IndexType left = 0;
1256       IndexType right = count - 1;
1257       for (;;)
1258       {
1259         while (left <= right && dataset_get(ind[left], cutfeat) < cutval) ++left;
1260         while (right && left <= right && dataset_get(ind[right], cutfeat) >= cutval) --right;
1261         if (left > right || !right) break;  // "!right" was added to support unsigned Index types
1262         std::swap(ind[left], ind[right]);
1263         ++left;
1264         --right;
1265       }
1266       /* If either list is empty, it means that all remaining features
1267              * are identical. Split in the middle to maintain a balanced tree.
1268              */
1269       lim1 = left;
1270       right = count - 1;
1271       for (;;)
1272       {
1273         while (left <= right && dataset_get(ind[left], cutfeat) <= cutval) ++left;
1274         while (right && left <= right && dataset_get(ind[right], cutfeat) > cutval) --right;
1275         if (left > right || !right) break;  // "!right" was added to support unsigned Index types
1276         std::swap(ind[left], ind[right]);
1277         ++left;
1278         --right;
1279       }
1280       lim2 = left;
1281     }
1282 
1283     DistanceType computeInitialDistances(const ElementType* vec, distance_vector_t& dists) const
1284     {
1285       assert(vec);
1286       DistanceType distsq = DistanceType();
1287 
1288       for (int i = 0; i < (DIM > 0 ? DIM : dim); ++i)
1289       {
1290         if (vec[i] < root_bbox[i].low)
1291         {
1292           dists[i] = distance.accum_dist(vec[i], root_bbox[i].low, i);
1293           distsq += dists[i];
1294         }
1295         if (vec[i] > root_bbox[i].high)
1296         {
1297           dists[i] = distance.accum_dist(vec[i], root_bbox[i].high, i);
1298           distsq += dists[i];
1299         }
1300       }
1301 
1302       return distsq;
1303     }
1304 
1305     /**
1306          * Performs an exact search in the tree starting from a node.
1307          * \tparam RESULTSET Should be any ResultSet<DistanceType>
1308          */
1309     template <class RESULTSET>
1310     void searchLevel(RESULTSET& result_set, const ElementType* vec, const NodePtr node, DistanceType mindistsq,
1311                      distance_vector_t& dists, const float epsError) const
1312     {
1313       /* If this is a leaf node, then do check and return. */
1314       if ((node->child1 == NULL) && (node->child2 == NULL))
1315       {
1316         //count_leaf += (node->lr.right-node->lr.left);  // Removed since was neither used nor returned to the user.
1317         DistanceType worst_dist = result_set.worstDist();
1318         for (IndexType i = node->node_type.lr.left; i < node->node_type.lr.right; ++i)
1319         {
1320           const IndexType index = vind[i];  // reorder... : i;
1321           DistanceType dist = distance(vec, index, (DIM > 0 ? DIM : dim));
1322           if (dist < worst_dist)
1323           {
1324             result_set.addPoint(dist, vind[i]);
1325           }
1326         }
1327         return;
1328       }
1329 
1330       /* Which child branch should be taken first? */
1331       int idx = node->node_type.sub.divfeat;
1332       ElementType val = vec[idx];
1333       DistanceType diff1 = val - node->node_type.sub.divlow;
1334       DistanceType diff2 = val - node->node_type.sub.divhigh;
1335 
1336       NodePtr bestChild;
1337       NodePtr otherChild;
1338       DistanceType cut_dist;
1339       if ((diff1 + diff2) < 0)
1340       {
1341         bestChild = node->child1;
1342         otherChild = node->child2;
1343         cut_dist = distance.accum_dist(val, node->node_type.sub.divhigh, idx);
1344       }
1345       else
1346       {
1347         bestChild = node->child2;
1348         otherChild = node->child1;
1349         cut_dist = distance.accum_dist(val, node->node_type.sub.divlow, idx);
1350       }
1351 
1352       /* Call recursively to search next level down. */
1353       searchLevel(result_set, vec, bestChild, mindistsq, dists, epsError);
1354 
1355       DistanceType dst = dists[idx];
1356       mindistsq = mindistsq + cut_dist - dst;
1357       dists[idx] = cut_dist;
1358       if (mindistsq * epsError <= result_set.worstDist())
1359       {
1360         searchLevel(result_set, vec, otherChild, mindistsq, dists, epsError);
1361       }
1362       dists[idx] = dst;
1363     }
1364 
1365    public:
1366     /**  Stores the index in a binary file.
1367           *   IMPORTANT NOTE: The set of data points is NOT stored in the file, so when loading the index object it must be constructed associated to the same source of data points used while building it.
1368           * See the example: examples/saveload_example.cpp
1369           * \sa loadIndex  */
1370     void saveIndex(FILE* stream)
1371     {
1372       save_value(stream, m_size);
1373       save_value(stream, dim);
1374       save_value(stream, root_bbox);
1375       save_value(stream, m_leaf_max_size);
1376       save_value(stream, vind);
1377       save_tree(stream, root_node);
1378     }
1379 
1380     /**  Loads a previous index from a binary file.
1381           *   IMPORTANT NOTE: The set of data points is NOT stored in the file, so the index object must be constructed associated to the same source of data points used while building the index.
1382           * See the example: examples/saveload_example.cpp
1383           * \sa loadIndex  */
1384     void loadIndex(FILE* stream)
1385     {
1386       load_value(stream, m_size);
1387       load_value(stream, dim);
1388       load_value(stream, root_bbox);
1389       load_value(stream, m_leaf_max_size);
1390       load_value(stream, vind);
1391       load_tree(stream, root_node);
1392     }
1393 
1394   };  // class KDTree
1395 
1396   /** An L2-metric KD-tree adaptor for working with data directly stored in an Eigen Matrix, without duplicating the data storage.
1397       *  Each row in the matrix represents a point in the state space.
1398       *
1399       *  Example of usage:
1400       * \code
1401       *     Eigen::Matrix<num_t,Dynamic,Dynamic>  mat;
1402       *     // Fill out "mat"...
1403       *
1404       *     typedef KDTreeEigenMatrixAdaptor< Eigen::Matrix<num_t,Dynamic,Dynamic> >  my_kd_tree_t;
1405       *     const int max_leaf = 10;
1406       *     my_kd_tree_t   mat_index(dimdim, mat, max_leaf );
1407       *     mat_index.index->buildIndex();
1408       *     mat_index.index->...
1409       * \endcode
1410       *
1411       *  \tparam DIM If set to >0, it specifies a compile-time fixed dimensionality for the points in the data set, allowing more compiler optimizations.
1412       *  \tparam Distance The distance metric to use: nanoflann::metric_L1, nanoflann::metric_L2, nanoflann::metric_L2_Simple, etc.
1413       */
1414   template <class MatrixType, int DIM = -1, class Distance = nanoflann::metric_L2>
1415   struct KDTreeEigenMatrixAdaptor
1416   {
1417     typedef KDTreeEigenMatrixAdaptor<MatrixType, DIM, Distance> self_t;
1418     typedef typename MatrixType::Scalar num_t;
1419     typedef typename MatrixType::Index IndexType;
1420     typedef typename Distance::template traits<num_t, self_t>::distance_t metric_t;
1421     typedef KDTreeSingleIndexAdaptor<metric_t, self_t, DIM, IndexType> index_t;
1422 
1423     index_t* index;  //! The kd-tree index for the user to call its methods as usual with any other FLANN index.
1424 
1425     /// Constructor: takes a const ref to the matrix object with the data points
1426     KDTreeEigenMatrixAdaptor(const int dimensionality, const MatrixType& mat, const int leaf_max_size = 10)
1427       : m_data_matrix(mat)
1428     {
1429       const IndexType dims = mat.cols();
1430       if (dims != dimensionality) throw std::runtime_error("Error: 'dimensionality' must match column count in data matrix");
1431       if (DIM > 0 && static_cast<int>(dims) != DIM)
1432         throw std::runtime_error("Data set dimensionality does not match the 'DIM' template argument");
1433       // cppcheck-suppress noOperatorEq
1434       index = new index_t(dims, *this /* adaptor */, nanoflann::KDTreeSingleIndexAdaptorParams(leaf_max_size));
1435       index->buildIndex();
1436     }
1437 
1438    private:
1439     /** Hidden copy constructor, to disallow copying this class (Not implemented) */
1440     KDTreeEigenMatrixAdaptor(const self_t&);
1441 
1442    public:
1443     ~KDTreeEigenMatrixAdaptor()
1444     {
1445       delete index;
1446     }
1447 
1448     const MatrixType& m_data_matrix;
1449 
1450     /** Query for the \a num_closest closest points to a given point (entered as query_point[0:dim-1]).
1451           *  Note that this is a short-cut method for index->findNeighbors().
1452           *  The user can also call index->... methods as desired.
1453           * \note nChecks_IGNORED is ignored but kept for compatibility with the original FLANN interface.
1454           */
1455     inline void query(const num_t* query_point, const size_t num_closest, IndexType* out_indices, num_t* out_distances_sq, const int /* nChecks_IGNORED */ = 10) const
1456     {
1457       nanoflann::KNNResultSet<num_t, IndexType> resultSet(num_closest);
1458       resultSet.init(out_indices, out_distances_sq);
1459       index->findNeighbors(resultSet, query_point, nanoflann::SearchParams());
1460     }
1461 
1462     /** @name Interface expected by KDTreeSingleIndexAdaptor
1463           * @{ */
1464 
1465     const self_t& derived() const
1466     {
1467       return *this;
1468     }
1469     self_t& derived()
1470     {
1471       return *this;
1472     }
1473 
1474     // Must return the number of data points
1475     inline size_t kdtree_get_point_count() const
1476     {
1477       return m_data_matrix.rows();
1478     }
1479 
1480     // Returns the L2 distance between the vector "p1[0:size-1]" and the data point with index "idx_p2" stored in the class:
1481     inline num_t kdtree_distance(const num_t* p1, const IndexType idx_p2, IndexType size) const
1482     {
1483       num_t s = 0;
1484       for (IndexType i = 0; i < size; i++)
1485       {
1486         const num_t d = p1[i] - m_data_matrix.coeff(idx_p2, i);
1487         s += d * d;
1488       }
1489       return s;
1490     }
1491 
1492     // Returns the dim'th component of the idx'th point in the class:
1493     inline num_t kdtree_get_pt(const IndexType idx, int dim) const
1494     {
1495       return m_data_matrix.coeff(idx, IndexType(dim));
1496     }
1497 
1498     // Optional bounding-box computation: return false to default to a standard bbox computation loop.
1499     //   Return true if the BBOX was already computed by the class and returned in "bb" so it can be avoided to redo it again.
1500     //   Look at bb.size() to find out the expected dimensionality (e.g. 2 or 3 for point clouds)
1501     template <class BBOX>
1502     bool kdtree_get_bbox(BBOX& /*bb*/) const
1503     {
1504       return false;
1505     }
1506 
1507     /** @} */
1508 
1509   };  // end of KDTreeEigenMatrixAdaptor
1510       /** @} */
1511 
1512   /** @} */  // end of grouping
1513 }  // namespace nanoflann
1514 
1515 #endif /* NANOFLANN_HPP_ */