Change neighborhood API, using a Concept + models

This commit is contained in:
Simon Giraudot 2016-11-03 12:03:45 +01:00
parent 31278c59cb
commit ed583ac6e9
7 changed files with 123 additions and 95 deletions

View File

@ -29,13 +29,16 @@
\cgalClassifedRefPages
## Concepts ##
- `CGAL::NeighborQuery`
## Classification ##
- `CGAL::Classifier<RandomAccessIterator, ItemMap>`
## Data Structures ##
- `CGAL::Classification::Helper<Kernel, RandomAccessIterator, PointMap, DiagonalizeTraits>`
- `CGAL::Classification::Neighborhood<Kernel, RandomAccessIterator, PointMap>`
- `CGAL::Classification::Point_set_neighborhood<Kernel, RandomAccessIterator, PointMap>`
- `CGAL::Classification::Local_eigen_analysis<Kernel, RandomAccessIterator, PointMap, DiagonalizeTraits>`
- `CGAL::Classification::Planimetric_grid<Kernel, RandomAccessIterator, PointMap>`

View File

@ -114,7 +114,7 @@ int main (int argc, char** argv)
std::cerr << "Training" << std::endl;
psc.training(800); // 800 trials
psc.run_with_graphcut (helper.neighborhood(), 0.5);
psc.run_with_graphcut (helper.neighborhood().k_neighbor_query(12), 0.5);
// Save the output in a colored PLY format
std::ofstream f ("classification.ply");

View File

@ -5,6 +5,7 @@
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Classifier.h>
#include <CGAL/Classification/Point_set_neighborhood.h>
#include <CGAL/Classification/Planimetric_grid.h>
#include <CGAL/Classification/Attribute.h>
#include <CGAL/Classification/Attributes_eigen.h>
@ -22,9 +23,9 @@ typedef CGAL::Identity_property_map<Point> Pmap;
typedef CGAL::Classifier<Iterator, Pmap> Classification;
typedef CGAL::Classification::Planimetric_grid<Kernel, Iterator, Pmap> Planimetric_grid;
typedef CGAL::Classification::Neighborhood<Kernel, Iterator, Pmap> Neighborhood;
typedef CGAL::Classification::Local_eigen_analysis<Kernel, Iterator, Pmap> Local_eigen_analysis;
typedef CGAL::Classification::Planimetric_grid<Kernel, Iterator, Pmap> Planimetric_grid;
typedef CGAL::Classification::Point_set_neighborhood<Kernel, Iterator, Pmap> Neighborhood;
typedef CGAL::Classification::Local_eigen_analysis<Kernel, Iterator, Pmap> Local_eigen_analysis;
typedef CGAL::Classification::Type_handle Type_handle;
typedef CGAL::Classification::Attribute_handle Attribute_handle;
@ -66,7 +67,7 @@ int main (int argc, char** argv)
Planimetric_grid grid (pts.begin(), pts.end(), Pmap(), bbox, grid_resolution);
Neighborhood neighborhood (pts.begin(), pts.end(), Pmap());
double garbage;
Local_eigen_analysis eigen (pts.begin(), pts.end(), Pmap(), neighborhood, 6, garbage);
Local_eigen_analysis eigen (pts.begin(), pts.end(), Pmap(), neighborhood.k_neighbor_query(6), garbage);
//! [Analysis]
///////////////////////////////////////////////////////////////////
@ -147,7 +148,7 @@ int main (int argc, char** argv)
///////////////////////////////////////////////////////////////////
// Run classification
psc.run_with_graphcut (neighborhood, 0.2);
psc.run_with_graphcut (neighborhood.k_neighbor_query(12), 0.2);
// Save the output in a colored PLY format

View File

@ -22,6 +22,7 @@
#include <CGAL/property_map.h>
#include <CGAL/Classifier.h>
#include <CGAL/Classification/Point_set_neighborhood.h>
#include <CGAL/Classification/Planimetric_grid.h>
#include <CGAL/Classification/Local_eigen_analysis.h>
#include <CGAL/Classification/Attribute.h>
@ -88,7 +89,7 @@ public:
<RandomAccessIterator, PointMap> Classifier;
typedef CGAL::Classification::Planimetric_grid
<Kernel, RandomAccessIterator, PointMap> Planimetric_grid;
typedef CGAL::Classification::Neighborhood
typedef CGAL::Classification::Point_set_neighborhood
<Kernel, RandomAccessIterator, PointMap> Neighborhood;
typedef CGAL::Classification::Local_eigen_analysis
<Kernel, RandomAccessIterator, PointMap, DiagonalizeTraits> Local_eigen_analysis;
@ -157,7 +158,7 @@ private:
t.reset();
t.start();
double range;
eigen = new Local_eigen_analysis (begin, end, point_map, *neighborhood, (std::size_t)6, range);
eigen = new Local_eigen_analysis (begin, end, point_map, neighborhood->k_neighbor_query(6), range);
if (this->voxel_size < 0)
this->voxel_size = range / 3;
t.stop();

View File

@ -9,8 +9,6 @@
#include <CGAL/Default_diagonalize_traits.h>
#include <CGAL/centroid.h>
#include <CGAL/Classification/Neighborhood.h>
namespace CGAL {
namespace Classification {
@ -37,7 +35,6 @@ public:
typedef typename Kernel::Point_3 Point;
typedef typename Kernel::Vector_3 Vector;
typedef typename Kernel::Plane_3 Plane;
typedef Classification::Neighborhood<Kernel, RandomAccessIterator, PointMap> Neighborhood;
/// \endcond
typedef CGAL::cpp11::array<double, 3> Eigenvalues; ///< Eigenvalues (sorted in ascending order)
@ -60,63 +57,25 @@ public:
/// \endcond
/*!
\brief Computes the local eigen analysis of an input range based
on a fixed number of local neighbors.
\brief Computes the local eigen analysis of an input range based
on a local neighborhood.
\tparam NeighborQuery is a model of `NeighborQuery`
\param begin Iterator to the first input object
\param end Past-the-end iterator
\param point_map Property map to access the input points
\param neighborhood Object used to access neighborhoods of points
\param radius_neighbors Radius of the local neighborhood
*/
Local_eigen_analysis (RandomAccessIterator begin,
RandomAccessIterator end,
PointMap point_map,
Neighborhood& neighborhood,
double radius_neighbors)
{
std::size_t size = end - begin;
m_eigenvalues.reserve (size);
m_centroids.reserve (size);
m_smallest_eigenvectors.reserve (size);
m_middle_eigenvectors.reserve (size);
m_largest_eigenvectors.reserve (size);
double nb = 0.;
for (std::size_t i = 0; i < size; i++)
{
std::vector<std::size_t> neighbors;
neighborhood.range_neighbors (get(point_map, begin[i]), radius_neighbors, std::back_inserter (neighbors));
nb += neighbors.size();
std::vector<Point> neighbor_points;
for (std::size_t j = 0; j < neighbors.size(); ++ j)
neighbor_points.push_back (get(point_map, begin[neighbors[j]]));
compute (get(point_map, begin[i]), neighbor_points);
}
std::cerr << "Mean number of nearest neighbors: " << nb / size << std::endl;
}
/*!
\brief Computes the local eigen analysis of an input range based
on a fixed radius local neighborhood.
\param begin Iterator to the first input object
\param end Past-the-end iterator
\param point_map Property map to access the input points
\param neighborhood Object used to access neighborhoods of points
\param knn Number of nearest neighbors used
\param neighbor_query Object used to access neighborhoods of points
\param mean_range The mean value of the range corresponding to the
`knn` number of neighbors is returned by the constructor through
this reference.
*/
template <typename NeighborQuery>
Local_eigen_analysis (RandomAccessIterator begin,
RandomAccessIterator end,
PointMap point_map,
Neighborhood& neighborhood,
std::size_t knn,
const NeighborQuery& neighbor_query,
double& mean_range)
{
std::size_t size = end - begin;
@ -131,7 +90,7 @@ public:
for (std::size_t i = 0; i < size; i++)
{
std::vector<std::size_t> neighbors;
neighborhood.k_neighbors (get(point_map, begin[i]), knn, std::back_inserter (neighbors));
neighbor_query (get(point_map, begin[i]), std::back_inserter (neighbors));
std::vector<Point> neighbor_points;
for (std::size_t j = 0; j < neighbors.size(); ++ j)

View File

@ -1,5 +1,5 @@
#ifndef CGAL_CLASSIFICATION_NEIGHBORHOOD_H
#define CGAL_CLASSIFICATION_NEIGHBORHOOD_H
#ifndef CGAL_CLASSIFICATION_POINT_SET_NEIGHBORHOOD_H
#define CGAL_CLASSIFICATION_POINT_SET_NEIGHBORHOOD_H
#include <vector>
@ -29,8 +29,9 @@ namespace Classification {
\tparam PointMap is a model of `ReadablePropertyMap` with value type `Point_3<Kernel>`.
*/
template <typename Kernel, typename RandomAccessIterator, typename PointMap>
class Neighborhood
class Point_set_neighborhood
{
typedef typename Kernel::FT FT;
typedef typename Kernel::Point_3 Point;
@ -67,8 +68,71 @@ class Neighborhood
public:
/// \cond SKIP_IN_MANUAL
Neighborhood () : m_tree (NULL) { }
/*!
Functor that returns a neighborhood with a fixed number of points.
\cgalModels NeighborQuery
*/
class K_neighbor_query
{
public:
typedef Point_set_neighborhood::Point value_type; ///<
private:
const Point_set_neighborhood& neighborhood;
std::size_t k;
public:
/*!
\brief Constructs a K neighbor query object.
\param neighborhood The point set neighborhood structure.
\param k The number of neighbors per query.
*/
K_neighbor_query (const Point_set_neighborhood& neighborhood, std::size_t k)
: neighborhood (neighborhood), k(k) { }
/// \cond SKIP_IN_MANUAL
template <typename OutputIterator>
void operator() (const value_type& query, OutputIterator output) const
{
neighborhood.k_neighbors (query, k, output);
}
/// \endcond
};
/*!
Functor that returns a neighborhood with a fixed radius.
\cgalModels NeighborQuery
*/
class Range_neighbor_query
{
public:
typedef Point_set_neighborhood::Point value_type; ///<
private:
const Point_set_neighborhood& neighborhood;
double radius;
public:
/*!
\brief Constructs a range neighbor query object.
\param neighborhood The point set neighborhood structure.
\param radius The radius of the neighbor query.
*/
Range_neighbor_query (const Point_set_neighborhood& neighborhood, double radius)
: neighborhood (neighborhood), radius(radius) { }
/// \cond SKIP_IN_MANUAL
template <typename OutputIterator>
void operator() (const value_type& query, OutputIterator output) const
{
neighborhood.range_neighbors (query, radius, output);
}
/// \endcond
};
friend class K_neighbor_query;
friend class Range_neighbor_query;
/// \cond SKIP_IN_MANUAL
Point_set_neighborhood () : m_tree (NULL) { }
/// \endcond
/*!
@ -78,7 +142,7 @@ public:
\param end Past-the-end iterator
\param point_map Property map to access the input points
*/
Neighborhood (const RandomAccessIterator& begin,
Point_set_neighborhood (const RandomAccessIterator& begin,
const RandomAccessIterator& end,
PointMap point_map)
: m_tree (NULL)
@ -107,7 +171,7 @@ public:
\param point_map Property map to access the input points
\param voxel_size
*/
Neighborhood (const RandomAccessIterator& begin,
Point_set_neighborhood (const RandomAccessIterator& begin,
const RandomAccessIterator& end,
PointMap point_map,
double voxel_size)
@ -130,7 +194,7 @@ public:
}
/// \cond SKIP_IN_MANUAL
~Neighborhood ()
~Point_set_neighborhood ()
{
if (m_tree != NULL)
delete m_tree;
@ -138,13 +202,23 @@ public:
/// \endcond
/*!
\brief Gets the nearest neighbors computed in a local sphere of user defined radius.
\tparam OutputIterator Iterator that must point to values of type `std::size_t`.
\param query The query point.
\param radius_neighbors Radius of the query sphere.
\param output Where the indices of found neighbor points are stored.
\brief Returns a neighbor query object with fixed number of neighbors `k`.
*/
K_neighbor_query k_neighbor_query (const std::size_t k) const
{
return K_neighbor_query (*this, k);
}
/*!
\brief Returns a neighbor query object with fixed radius `radius`.
*/
Range_neighbor_query range_neighbor_query (const double radius) const
{
return Range_neighbor_query (*this, radius);
}
private:
template <typename OutputIterator>
void range_neighbors (const Point& query, const FT radius_neighbors, OutputIterator output) const
{
@ -153,14 +227,6 @@ public:
m_tree->search (output, fs);
}
/*!
\brief Gets the K nearest neighbors.
\tparam OutputIterator Iterator that must point to values of type `std::size_t`.
\param query The query point.
\param k Number of nearest neighbors.
\param output Where the indices of found neighbor points are stored.
*/
template <typename OutputIterator>
void k_neighbors (const Point& query, const std::size_t k, OutputIterator output) const
{
@ -170,7 +236,6 @@ public:
*(output ++) = it->first;
}
private:
/// \cond SKIP_IN_MANUAL
template <typename Map>
void voxelize_point_set (std::vector<std::size_t>& indices, Map point_map,
@ -223,4 +288,4 @@ private:
}
#endif // CGAL_CLASSIFICATION_NEIGHBORHOOD_H
#endif // CGAL_CLASSIFICATION_POINT_SET_POINT_SET_NEIGHBORHOOD_H

View File

@ -231,12 +231,11 @@ public:
local neighborhood of items. This method is a compromise between
efficiency and reliability.
\param neighborhood Object used to access neighborhoods of items
\param radius_neighbors Radius used for smoothing
\tparam NeighborQuery is a model of `NeighborQuery`
\param neighbor_query is used to access neighborhoods of items
*/
template <typename Neighborhood>
void run_with_local_smoothing (const Neighborhood& neighborhood,
const double radius_neighbors)
template <typename NeighborQuery>
void run_with_local_smoothing (const NeighborQuery& neighbor_query)
{
prepare_classification ();
@ -250,8 +249,7 @@ public:
for (std::size_t s=0; s < m_input.size(); ++ s)
{
std::vector<std::size_t> neighbors;
neighborhood.range_neighbors (m_input[s], radius_neighbors,
std::back_inserter (neighbors));
neighbor_query (m_input[s], std::back_inserter (neighbors));
std::vector<double> mean (values.size(), 0.);
for (std::size_t n = 0; n < neighbors.size(); ++ n)
@ -296,15 +294,16 @@ public:
and alpha-expansion algorithm. This method is slow but provides
the user with good quality results.
\param neighborhood Object used to access neighborhoods of items
\tparam NeighborQuery is a model of `NeighborQuery`
\param neighbor_query is used to access neighborhoods of items
\param weight Weight of the regularization with respect to the
classification energy. Higher values produce more regularized
output but may result in a loss of details.
*/
template <typename Neighborhood>
void run_with_graphcut (const Neighborhood& neighborhood,
const double weight = 0.5)
template <typename NeighborQuery>
void run_with_graphcut (const NeighborQuery& neighbor_query,
const double& weight)
{
prepare_classification ();
@ -321,7 +320,7 @@ public:
{
std::vector<std::size_t> neighbors;
neighborhood.k_neighbors (m_input[s], 12, std::back_inserter (neighbors));
neighbor_query (m_input[s], std::back_inserter (neighbors));
for (std::size_t i = 0; i < neighbors.size(); ++ i)
if (s != neighbors[i])