Skip to content

Commit

Permalink
Use std::cref() for passing Eigen matrix datasets
Browse files Browse the repository at this point in the history
  • Loading branch information
jlblancoc committed Jul 27, 2018
1 parent df61b0b commit ad7547f
Show file tree
Hide file tree
Showing 3 changed files with 87 additions and 75 deletions.
1 change: 1 addition & 0 deletions CHANGELOG.txt
Original file line number Diff line number Diff line change
@@ -1,4 +1,5 @@
nanoflann 1.3.0: (Unreleased)
* Eigen::Matrix datasets: now uses std::cref() to store a reference to matrix.
* GSOC2017 contributions by Pranjal Kumar Rai:
* Support for dynamic datasets.
* Support for non-Euclidean spaces: SO(2), SO(3)
Expand Down
144 changes: 76 additions & 68 deletions examples/matrix_example.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -28,8 +28,8 @@

#include <nanoflann.hpp>

#include <ctime>
#include <cstdlib>
#include <ctime>
#include <iostream>

#include <Eigen/Dense>
Expand All @@ -41,77 +41,85 @@ using namespace nanoflann;
const int SAMPLES_DIM = 15;

template <typename Der>
void generateRandomPointCloud(Eigen::MatrixBase<Der> &mat, const size_t N, const size_t dim, const typename Der::Scalar max_range = 10)
{
std::cout << "Generating "<< N << " random points...";
mat.resize(N,dim);
for (size_t i = 0; i < N; i++)
for (size_t d = 0; d < dim; d++)
mat(i,d) = max_range * (rand() % 1000) / typename Der::Scalar(1000);
std::cout << "done\n";
void generateRandomPointCloud(Eigen::MatrixBase<Der> &mat, const size_t N,
const size_t dim,
const typename Der::Scalar max_range = 10) {
std::cout << "Generating " << N << " random points...";
mat.resize(N, dim);
for (size_t i = 0; i < N; i++)
for (size_t d = 0; d < dim; d++)
mat(i, d) = max_range * (rand() % 1000) / typename Der::Scalar(1000);
std::cout << "done\n";
}

template <typename num_t>
void kdtree_demo(const size_t nSamples, const size_t dim)
{
Eigen::Matrix<num_t, Dynamic, Dynamic> mat(nSamples, dim);

const num_t max_range = 20;

// Generate points:
generateRandomPointCloud(mat, nSamples, dim, max_range);

// cout << mat << endl;

// Query point:
std::vector<num_t> query_pt(dim);
for (size_t d = 0; d < dim; d++)
query_pt[d] = max_range * (rand() % 1000) / num_t(1000);


// ------------------------------------------------------------
// construct a kd-tree index:
// Some of the different possibilities (uncomment just one)
// ------------------------------------------------------------
// Dimensionality set at run-time (default: L2)
typedef KDTreeEigenMatrixAdaptor< Eigen::Matrix<num_t, Dynamic, Dynamic> > my_kd_tree_t;

// Dimensionality set at compile-time
// typedef KDTreeEigenMatrixAdaptor< Eigen::Matrix<num_t,Dynamic,Dynamic> > my_kd_tree_t;

// Dimensionality set at compile-time: Explicit selection of the distance metric: L2
// typedef KDTreeEigenMatrixAdaptor< Eigen::Matrix<num_t,Dynamic,Dynamic>,nanoflann::metric_L2> my_kd_tree_t;

// Dimensionality set at compile-time: Explicit selection of the distance metric: L2_simple
// typedef KDTreeEigenMatrixAdaptor< Eigen::Matrix<num_t,Dynamic,Dynamic>,nanoflann::metric_L2_Simple> my_kd_tree_t;

// Dimensionality set at compile-time: Explicit selection of the distance metric: L1
// typedef KDTreeEigenMatrixAdaptor< Eigen::Matrix<num_t,Dynamic,Dynamic>,nanoflann::metric_L1> my_kd_tree_t;

my_kd_tree_t mat_index(mat, 10 /* max leaf */ );
mat_index.index->buildIndex();

// do a knn search
const size_t num_results = 3;
vector<size_t> ret_indexes(num_results);
vector<num_t> out_dists_sqr(num_results);

nanoflann::KNNResultSet<num_t> resultSet(num_results);

resultSet.init(&ret_indexes[0], &out_dists_sqr[0] );
mat_index.index->findNeighbors(resultSet, &query_pt[0], nanoflann::SearchParams(10));

std::cout << "knnSearch(nn="<<num_results<<"): \n";
for (size_t i=0; i<num_results; i++)
std::cout << "ret_index["<<i<<"]=" << ret_indexes[i] << " out_dist_sqr=" << out_dists_sqr[i] << endl;

void kdtree_demo(const size_t nSamples, const size_t dim) {
Eigen::Matrix<num_t, Dynamic, Dynamic> mat(nSamples, dim);

const num_t max_range = 20;

// Generate points:
generateRandomPointCloud(mat, nSamples, dim, max_range);

// cout << mat << endl;

// Query point:
std::vector<num_t> query_pt(dim);
for (size_t d = 0; d < dim; d++)
query_pt[d] = max_range * (rand() % 1000) / num_t(1000);

// ------------------------------------------------------------
// construct a kd-tree index:
// Some of the different possibilities (uncomment just one)
// ------------------------------------------------------------
// Dimensionality set at run-time (default: L2)
typedef KDTreeEigenMatrixAdaptor<Eigen::Matrix<num_t, Dynamic, Dynamic>>
my_kd_tree_t;

// Dimensionality set at compile-time
// typedef KDTreeEigenMatrixAdaptor< Eigen::Matrix<num_t,Dynamic,Dynamic> >
// my_kd_tree_t;

// Dimensionality set at compile-time: Explicit selection of the distance
// metric: L2
// typedef KDTreeEigenMatrixAdaptor<
// Eigen::Matrix<num_t,Dynamic,Dynamic>,nanoflann::metric_L2> my_kd_tree_t;

// Dimensionality set at compile-time: Explicit selection of the distance
// metric: L2_simple
// typedef KDTreeEigenMatrixAdaptor<
// Eigen::Matrix<num_t,Dynamic,Dynamic>,nanoflann::metric_L2_Simple>
// my_kd_tree_t;

// Dimensionality set at compile-time: Explicit selection of the distance
// metric: L1
// typedef KDTreeEigenMatrixAdaptor<
// Eigen::Matrix<num_t,Dynamic,Dynamic>,nanoflann::metric_L1> my_kd_tree_t;

my_kd_tree_t mat_index(dim, std::cref(mat), 10 /* max leaf */);
mat_index.index->buildIndex();

// do a knn search
const size_t num_results = 3;
vector<size_t> ret_indexes(num_results);
vector<num_t> out_dists_sqr(num_results);

nanoflann::KNNResultSet<num_t> resultSet(num_results);

resultSet.init(&ret_indexes[0], &out_dists_sqr[0]);
mat_index.index->findNeighbors(resultSet, &query_pt[0],
nanoflann::SearchParams(10));

std::cout << "knnSearch(nn=" << num_results << "): \n";
for (size_t i = 0; i < num_results; i++)
std::cout << "ret_index[" << i << "]=" << ret_indexes[i]
<< " out_dist_sqr=" << out_dists_sqr[i] << endl;
}

int main(int argc, char** argv)
{
// Randomize Seed
srand(time(NULL));
kdtree_demo<float>(1e3 /* samples */, SAMPLES_DIM /* dim */);
int main(int argc, char **argv) {
// Randomize Seed
srand(time(NULL));
kdtree_demo<float>(1e3 /* samples */, SAMPLES_DIM /* dim */);

return 0;
return 0;
}
17 changes: 10 additions & 7 deletions include/nanoflann.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -1704,10 +1704,11 @@ class KDTreeSingleIndexDynamicAdaptor_
if (dist < worst_dist) {
if (!result_set.addPoint(
static_cast<typename RESULTSET::DistanceType>(dist),
static_cast<typename RESULTSET::IndexType>(vind[i]))) {
static_cast<typename RESULTSET::IndexType>(
BaseClassRef::vind[i]))) {
// the resultset doesn't want to receive any more points, we're done
// searching!
return false;
return; // false;
}
}
}
Expand Down Expand Up @@ -1941,9 +1942,9 @@ class KDTreeSingleIndexDynamicAdaptor {
* Distance The distance metric to use: nanoflann::metric_L1,
* nanoflann::metric_L2, nanoflann::metric_L2_Simple, etc.
*/
template <class MatrixType, class Distance = nanoflann::metric_L2>
template <class MatrixType, int DIM = -1, class Distance = nanoflann::metric_L2>
struct KDTreeEigenMatrixAdaptor {
typedef KDTreeEigenMatrixAdaptor<MatrixType, Distance> self_t;
typedef KDTreeEigenMatrixAdaptor<MatrixType, DIM, Distance> self_t;
typedef typename MatrixType::Scalar num_t;
typedef typename MatrixType::Index IndexType;
typedef
Expand All @@ -1960,7 +1961,7 @@ struct KDTreeEigenMatrixAdaptor {
const std::reference_wrapper<const MatrixType> &mat,
const int leaf_max_size = 10)
: m_data_matrix(mat) {
const int dims = static_cast<int>(mat.cols());
const int dims = static_cast<int>(mat.get().cols());
if (dims != dimensionality)
throw std::runtime_error(
"Error: 'dimensionality' must match column count in data matrix");
Expand Down Expand Up @@ -2002,11 +2003,13 @@ struct KDTreeEigenMatrixAdaptor {
self_t &derived() { return *this; }

// Must return the number of data points
inline size_t kdtree_get_point_count() const { return m_data_matrix.rows(); }
inline size_t kdtree_get_point_count() const {
return m_data_matrix.get().rows();
}

// Returns the dim'th component of the idx'th point in the class:
inline num_t kdtree_get_pt(const IndexType idx, int dim) const {
return m_data_matrix.coeff(idx, IndexType(dim));
return m_data_matrix.get().coeff(idx, IndexType(dim));
}

// Optional bounding-box computation: return false to default to a standard
Expand Down

0 comments on commit ad7547f

Please sign in to comment.