\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 4.12 - Point Set Shape Detection
User Manual

Authors
Sven Oesau, Yannick Verdie, Clément Jamin, Pierre Alliez, Florent Lafarge, Simon Giraudot

Introduction

This CGAL component implements two algorithms for shape detection:

  • the efficient RANSAC method for shape detection, contributed by Schnabel et al. [2];
  • the region growing method, contributed by Lafarge and Mallet [1].

From an unstructured point set with unoriented normals, these algorithms detect a set of shapes (see Figure 72.1). Five types of primitive shapes are provided by this package: plane, sphere, cylinder, cone and torus. Other types of primitive shapes can also be added by the user.

Note that at the moment, region growing only detects planes.

overview2.png
Figure 72.1 (a) Input point set. (b) Point set depicted with one color per detected shape.

Methods

Both methods take as input a point set with unoriented normals and provide as output a set of detected shapes with associated input points. The output of the shape detection algorithms is a set of detected shapes with assigned points and all remaining points not covered by these shapes. Each input point can be assigned to at most one detected shape.

Efficient RANSAC

The shapes are detected via a RANSAC-type approach, i.e., a random sample consensus. The basic RANSAC approach repeats the following steps: (1) Randomly select samples from the input points. (2) Fit a shape to the selected samples. (3) Count the number of inliers to the shape, inliers being within a user-specified error tolerance to the shape. Steps 1-3 are repeated for a prescribed number of iterations and the shape with highest number of inliers, referred to as largest shape, is kept.

In our context, the error between a point and a shape is defined by its distance and normal deviation to the shape. A random subset corresponds to the minimal number of points (with normals) required to uniquely define a primitive.

For very large point sets the basic RANSAC method is not practical when testing all possible shape candidates against the input data in order to find the largest shape. The main idea behind the efficient RANSAC method consists in testing shape candidates against subsets of the input data. Shape candidates are constructed until the probability to miss the largest candidate is lower than a user-specified parameter. The largest shape is repeatedly extracted until no more shapes, restricted to cover a minimum number of points, can be extracted. An additional gain in efficiency is achieved through exploiting the normal attributes during initial shape construction and enumeration of inliers.

The support of a shape refers to the footprint of the points covered by the primitive. To avoid generating shapes with fragmented support we enforce a connectivity constraint by considering only one connected component, referred to as cluster, selected as the one covering the largest number of inliers. See Section Parameters for more details.

Region Growing

Planes are detected by growing planar regions from estimated planar seeds. First, points are sorted by a local planarity measure (the fitting quality returned by linear_least_squares_fitting_3()), such that points whose local neighborhood is the most planar are treated first.

The following steps are repeated: (1) Pick the next available point in the sorted list. (2) Initialize a plane tangent to this point and perpendicular to the normal vector of this point. (3) Compute the local neighborhood of the point with a user-specified radius. (4) Find inliers among these points, i.e. points that are within a user-specified error tolerance to the plane (and with a normal deviation within a user-specified tolerance to the normal of the plane). Steps 3-4 are repeated until no other inlier can be added. The plane is periodically reestimated by principal component analysis using all the inliers. If, after the growing, the plane has less inliers than a user-specified minimum number of points, it is discarded and the points are made available for other shapes to take once again.

Comparisons

The efficient RANSAC is a very quick method. Because it is based on randomness, it is not deterministic and some small shapes might be missed in the detection process.

Region growing usually takes longer than the efficient RANSAC but may provide better quality output in the presence of large scenes with numerous small details: as it iterates throughout all points of the scene, there are less chances of missing a shape. In addition, it is deterministic (for a given input and a given set of parameters, it always returns the same output, whereas as RANSAC uses randomness to detect shapes, the output of a RANSAC algorithm may vary). See figure Figure 72.2.

comparison.png
Figure 72.2 Comparison of efficient RANSAC and region growing. Top: input point set. Bottom left: output of efficient RANSAC, \(78\%\) of the shapes correctly detected in 8 seconds. Bottom right: output of region growing, \(100\%\) of the shapes detected in 15 seconds. Unassigned points are in black in both output images.

Parameters

The algorithms have four parameters that are common. They have the same effects on the output for both RANSAC and region growing.

  • epsilon and normal_threshold: The error between a point-with-normal \(p\) and a shape \(S\) is defined by its Euclidean distance and normal deviation to \(S\). The normal deviation is computed between the normal at \(p\) and the normal of \(S\) at the closest projection of \(p\) onto \(S\). Parameter epsilon defines the absolute maximum tolerance Euclidean distance between a point and a shape. A high value for epsilon leads to the detection of fewer large shapes and hence a less detailed detection. A low value for epsilon yields a more detailed detection, but may lead to either lower coverage or over-segmentation. Over-segmentation translates into detection of fragmented shapes when epsilon is within or below the noise level. When the input point set is made of free-form parts, a higher tolerance epsilon allows for detecting more primitive shapes that approximate some of the free-form surfaces. The impact of this parameter is depicted by Figure Figure 72.3. Its impact on performance is evaluated in Section Performances.
    epsilon_variation2.png
    Figure 72.3 Impact of epsilon over level of details of the detection. (a) Input point set. (b) Detection of planar shapes with epsilon set to 2.0 (one color per detected shape). Most details such as chimneys on the roof are not distinguished. (c) Detection with epsilon set to 0.5. The facades are correctly detected and some details of the roof are detected. (d) Setting epsilon to 0.25 yields a more detailed but slightly over-segmented detection.

  • cluster_epsilon: The region growing algorithm uses cluster_epsilon to compute the neighborhood around points when growing the regions while the efficient RANSAC uses this parameter to cluster the points into connected components covered by a detected shape. For developable shapes that admit a trivial planar parameterization (plane, cylinder, cone) the points covered by a shape are mapped to a 2D parameter space chosen to minimize distortion and best preserve arc-length distances. This 2D parameter space is discretized using a regular grid, and a connected component search is performed to identify the largest cluster. Parameter cluster_epsilon defines the spacing between two cells of the regular grid, so that two points separated by a distance of at most \(2\sqrt{2}\), cluster_epsilon are considered adjacent. For non-developable shapes the connected components are identified by computing a neighboring graph in 3D and walking in the graph. The impact of parameter cluster_epsilon is depicted by Figure Figure 72.4.

varying_connectivity.png
Figure 72.4 For efficient RANSAC, parameter cluster_epsilon controls the connectivity of the points covered by a detected shape. The input point set is sampled on four coplanar squares. (a) A large value for cluster_epsilon leads to detecting a single planar shape. (b) A moderate value for cluster_epsilon yields the detection of four squares. Notice that few points within the squares are not detected as not connected. (c) A small value for cluster_epsilon leads to over-segmentation.

  • min_points: This minimum number of points controls the termination of the algorithm. The shape search is iterated until no further shapes can be found with a higher support. Note that for efficient RANSAC, this parameter is not strict: depending on the chosen probability, shapes may be extracted with a number of points lower than the specified parameter.

In addition, the efficient RANSAC has a fifth parameter:

  • probability: This parameter defines the probability to miss the largest candidate shape. A lower probability provides a higher reliability and determinism at the cost of longer running times due to higher search endurance.

Examples

The main classes Shape_detection_3::Efficient_RANSAC and Shape_detection_3::Region_growing take a template parameter Shape_detection_3::Shape_detection_traits that defines the geometric types and input format. Property maps provide a means to interface with user-specific data structures. The first parameter of the Shape_detection_3::Shape_detection_traits class is the common Kernel. In order to match the constraints of property maps, an iterator type and two maps that map an iterator to a point and a normal are specified in the Shape_detection_3::Shape_detection_traits class. The concept behind property maps is detailed in CGAL and Propery Maps.

Typical usage consists of five steps:

  1. Provide input data via a range iterator
  2. Register shape factories
  3. Choose parameters
  4. Detect
  5. Retrieve detected shapes

Basic Planar Shape Detection

The following minimal example reads a point set from a file and detects only planar shapes. The default parameters are used for detection.


File Point_set_shape_detection_3/shape_detection_basic.cpp

#if defined (_MSC_VER) && !defined (_WIN64)
#pragma warning(disable:4244) // boost::number_distance::distance()
// converts 64 to 32 bits integers
#endif
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Point_with_normal_3.h>
#include <CGAL/property_map.h>
#include <iostream>
#include <fstream>
// Type declarations
typedef std::pair<Kernel::Point_3, Kernel::Vector_3> Point_with_normal;
typedef std::vector<Point_with_normal> Pwn_vector;
// In Shape_detection_traits the basic types, i.e., Point and Vector types
// as well as iterator type and property maps, are defined.
<Kernel, Pwn_vector, Point_map, Normal_map> Traits;
// This program both works for RANSAC and Region Growing
template <typename ShapeDetection>
int run(const char* filename)
{
// Points with normals.
Pwn_vector points;
// Loads point set from a file.
// read_xyz_points_and_normals takes an OutputIterator for storing the points
// and a property map to store the normal vector with each point.
std::ifstream stream(filename);
if (!stream ||
std::back_inserter(points),
CGAL::parameters::point_map(Point_map()).
normal_map(Normal_map())))
{
std::cerr << "Error: cannot read file cube.pwn" << std::endl;
return EXIT_FAILURE;
}
// Instantiates shape detection engine.
ShapeDetection shape_detection;
// Provides the input data.
shape_detection.set_input(points);
// Registers planar shapes via template method.
shape_detection.template add_shape_factory<Plane>();
// Detects registered shapes with default parameters.
shape_detection.detect();
// Prints number of detected shapes.
std::cout << shape_detection.shapes().end() - shape_detection.shapes().begin() << " shapes detected." << std::endl;
return EXIT_SUCCESS;
}
int main (int argc, char** argv)
{
if (argc > 1 && std::string(argv[1]) == "-r")
{
std::cout << "Efficient RANSAC" << std::endl;
return run<Efficient_ransac> ((argc > 2) ? argv[2] : "data/cube.pwn");
}
std::cout << "Region Growing" << std::endl;
return run<Region_growing> ((argc > 1) ? argv[1] : "data/cube.pwn");
}

Planar Shape Detection with Callback

Both Region Growing and Efficient RANSAC provide a callback mechanism that enables the user to track the progress of the algorithms. It can be used, for example, to terminate the algorithm based on a timeout. In the following example, the algorithm stops if it takes more than half a second and prints out the progress made:


File Point_set_shape_detection_3/shape_detection_with_callback.cpp

#if defined (_MSC_VER) && !defined (_WIN64)
#pragma warning(disable:4244) // boost::number_distance::distance()
// converts 64 to 32 bits integers
#endif
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Point_with_normal_3.h>
#include <CGAL/property_map.h>
#include <CGAL/Timer.h>
#include <iostream>
#include <fstream>
// Type declarations
typedef std::pair<Kernel::Point_3, Kernel::Vector_3> Point_with_normal;
typedef std::vector<Point_with_normal> Pwn_vector;
// In Shape_detection_traits the basic types, i.e., Point and Vector types
// as well as iterator type and property maps, are defined.
<Kernel, Pwn_vector, Point_map, Normal_map> Traits;
struct Timeout_callback
{
mutable int nb;
mutable CGAL::Timer timer;
const double limit;
Timeout_callback(double limit) : nb(0), limit(limit)
{
timer.start();
}
bool operator()(double advancement) const
{
// Avoid calling time() at every single iteration, which could
// impact performances very badly
++ nb;
if (nb % 1000 != 0)
return true;
// If the limit is reach, interrupt the algorithm
if (timer.time() > limit)
{
std::cerr << "Algorithm takes too long, exiting ("
<< 100. * advancement << "% done)" << std::endl;
return false;
}
return true;
}
};
// This program both works for RANSAC and Region Growing
template <typename ShapeDetection>
int run(const char* filename)
{
Pwn_vector points;
std::ifstream stream(filename);
if (!stream ||
std::back_inserter(points),
CGAL::parameters::point_map(Point_map()).
normal_map(Normal_map())))
{
std::cerr << "Error: cannot read file cube.pwn" << std::endl;
return EXIT_FAILURE;
}
ShapeDetection shape_detection;
shape_detection.set_input(points);
shape_detection.template add_shape_factory<Plane>();
// Create callback that interrupts the algorithm if it takes more than half a second
Timeout_callback timeout_callback(0.5);
// Detects registered shapes with default parameters.
shape_detection.detect(typename ShapeDetection::Parameters(),
timeout_callback);
return EXIT_SUCCESS;
}
int main (int argc, char** argv)
{
if (argc > 1 && std::string(argv[1]) == "-r")
{
std::cout << "Efficient RANSAC" << std::endl;
return run<Efficient_ransac> ((argc > 2) ? argv[2] : "data/cube.pwn");
}
std::cout << "Region Growing" << std::endl;
return run<Region_growing> ((argc > 1) ? argv[1] : "data/cube.pwn");
}

Setting Parameters and Using Different Types of Shape

This example illustrates the user selection of parameters using the Shape_detection_3::Efficient_RANSAC::Parameters class. Shape detection is performed on five types of shapes (plane, cylinder, sphere, cone, torus). The input point set is sampled on a surface mostly composed of piece-wise planar and cylindrical parts, in addition to free-form parts.

Basic information of the detected shapes is written to the standard output: if the plane is either a plane or a cylinder, specific parameters are recovered, otherwise the general method info() is used to get the shape parameters in a string object. Note that specific parameters can be recovered for any of the provided shape.


File Point_set_shape_detection_3/efficient_RANSAC_parameters.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Point_with_normal_3.h>
#include <CGAL/property_map.h>
#include <iostream>
#include <fstream>
// Type declarations
typedef Kernel::FT FT;
typedef std::pair<Kernel::Point_3, Kernel::Vector_3> Point_with_normal;
typedef std::vector<Point_with_normal> Pwn_vector;
// In Shape_detection_traits the basic types, i.e., Point and Vector types
// as well as iterator type and property maps, are defined.
Pwn_vector, Point_map, Normal_map> Traits;
int main()
{
// Points with normals.
Pwn_vector points;
// Loads point set from a file.
// read_xyz_points_and_normals takes an OutputIterator for storing the points
// and a property map to store the normal vector with each point.
std::ifstream stream("data/cube.pwn");
if (!stream ||
std::back_inserter(points),
CGAL::parameters::point_map(Point_map()).
normal_map(Normal_map())))
{
std::cerr << "Error: cannot read file cube.pwn" << std::endl;
return EXIT_FAILURE;
}
std::cout << points.size() << " points" << std::endl;
// Instantiates shape detection engine.
Efficient_ransac ransac;
// Provides the input data.
ransac.set_input(points);
// Register shapes for detection
ransac.add_shape_factory<Plane>();
ransac.add_shape_factory<Sphere>();
ransac.add_shape_factory<Cylinder>();
ransac.add_shape_factory<Cone>();
ransac.add_shape_factory<Torus>();
// Sets parameters for shape detection.
Efficient_ransac::Parameters parameters;
// Sets probability to miss the largest primitive at each iteration.
parameters.probability = 0.05;
// Detect shapes with at least 500 points.
parameters.min_points = 200;
// Sets maximum Euclidean distance between a point and a shape.
parameters.epsilon = 0.002;
// Sets maximum Euclidean distance between points to be clustered.
parameters.cluster_epsilon = 0.01;
// Sets maximum normal deviation.
// 0.9 < dot(surface_normal, point_normal);
parameters.normal_threshold = 0.9;
// Detects shapes
ransac.detect(parameters);
// Prints number of detected shapes and unassigned points.
std::cout << ransac.shapes().end() - ransac.shapes().begin() << " detected shapes, "
<< ransac.number_of_unassigned_points()
<< " unassigned points." << std::endl;
// Efficient_ransac::shapes() provides
// an iterator range to the detected shapes.
Efficient_ransac::Shape_range shapes = ransac.shapes();
Efficient_ransac::Shape_range::iterator it = shapes.begin();
while (it != shapes.end()) {
// Get specific parameters depending on detected shape.
if (Plane* plane = dynamic_cast<Plane*>(it->get()))
{
Kernel::Vector_3 normal = plane->plane_normal();
std::cout << "Plane with normal " << normal
<< std::endl;
// Plane shape can also be converted to Kernel::Plane_3
std::cout << "Kernel::Plane_3: " << static_cast<Kernel::Plane_3>(*plane) << std::endl;
}
else if (Cylinder* cyl = dynamic_cast<Cylinder*>(it->get()))
{
Kernel::Line_3 axis = cyl->axis();
FT radius = cyl->radius();
std::cout << "Cylinder with axis " << axis
<< " and radius " << radius
<< std::endl;
}
else
{
// Prints the parameters of the detected shape.
// This function is available for any type of shape.
std::cout << (*it)->info() << std::endl;
}
// Proceeds with next detected shape.
it++;
}
return EXIT_SUCCESS;
}

Retrieving Points Assigned to Shapes

This example illustrates how to access the points assigned to each shape and compute the mean error. A timer measures the running performance.


File Point_set_shape_detection_3/efficient_RANSAC_point_access.cpp

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Point_with_normal_3.h>
#include <CGAL/property_map.h>
#include <CGAL/Timer.h>
#include <CGAL/number_utils.h>
#include <iostream>
#include <fstream>
// Type declarations
typedef Kernel::FT FT;
typedef std::pair<Kernel::Point_3, Kernel::Vector_3> Point_with_normal;
typedef std::vector<Point_with_normal> Pwn_vector;
// In Shape_detection_traits the basic types, i.e., Point and Vector types
// as well as iterator type and property maps, are defined.
Pwn_vector, Point_map, Normal_map> Traits;
int main()
{
// Points with normals.
Pwn_vector points;
// Loads point set from a file.
// read_xyz_points_and_normals takes an OutputIterator for storing the points
// and a property map to store the normal vector with each point.
std::ifstream stream("data/cube.pwn");
if (!stream ||
std::back_inserter(points),
CGAL::parameters::point_map(Point_map()).
normal_map(Normal_map())))
{
std::cerr << "Error: cannot read file cube.pwn" << std::endl;
return EXIT_FAILURE;
}
// Instantiates shape detection engine.
Efficient_ransac ransac;
// Provides the input data.
ransac.set_input(points);
// Registers detection of planes
ransac.add_shape_factory<Plane>();
// Measures time before setting up the shape detection.
CGAL::Timer time;
time.start();
// Build internal data structures.
ransac.preprocess();
// Measures time after preprocessing.
time.stop();
std::cout << "preprocessing took: " << time.time() * 1000 << "ms" << std::endl;
// Perform detection several times and choose result with highest coverage.
Efficient_ransac::Shape_range shapes = ransac.shapes();
FT best_coverage = 0;
for (size_t i = 0;i<3;i++) {
// Reset timer.
time.reset();
time.start();
// Detects shapes.
ransac.detect();
// Measures time after detection.
time.stop();
// Compute coverage, i.e. ratio of the points assigned to a shape.
FT coverage = FT(points.size() - ransac.number_of_unassigned_points())
/ FT(points.size());
// Prints number of assigned shapes and unsassigned points.
std::cout << "time: " << time.time() * 1000 << "ms" << std::endl;
std::cout << ransac.shapes().end() - ransac.shapes().begin() << " primitives, "
<< coverage << " coverage" << std::endl;
// Choose result with highest coverage.
if (coverage > best_coverage) {
best_coverage = coverage;
// Efficient_ransac::shapes() provides
// an iterator range to the detected shapes.
shapes = ransac.shapes();
}
}
Efficient_ransac::Shape_range::iterator it = shapes.begin();
while (it != shapes.end()) {
boost::shared_ptr<Efficient_ransac::Shape> shape = *it;
// Using Shape_base::info() for printing
// the parameters of the detected shape.
std::cout << (*it)->info();
// Sums distances of points to detected shapes.
FT sum_distances = 0;
// Iterates through point indices assigned to each detected shape.
std::vector<std::size_t>::const_iterator
index_it = (*it)->indices_of_assigned_points().begin();
while (index_it != (*it)->indices_of_assigned_points().end()) {
// Retrieves point
const Point_with_normal &p = *(points.begin() + (*index_it));
// Adds Euclidean distance between point and shape.
sum_distances += CGAL::sqrt((*it)->squared_distance(p.first));
// Proceeds with next point.
index_it++;
}
// Computes and prints average distance.
FT average_distance = sum_distances / shape->indices_of_assigned_points().size();
std::cout << " average distance: " << average_distance << std::endl;
// Proceeds with next detected shape.
it++;
}
return EXIT_SUCCESS;
}

Adding More Shapes to RANSAC

Other types of shapes can be detected by implementing a shape class derived from the class Shape_base and registering it to the shape detection factory of the efficient RANSAC object. This class must provide the following functions: construct a shape from a small set of given points, compute the squared distance from a query point and the shape and compute the normal deviation between a query point with normal and the normal to the shape at the closest point from the query. The used shape parameters are added as members to the derived class.

Note that the RANSAC approach is efficient for shapes that are uniquely defined by a small number of points, denoted by the number of required samples. The algorithm aims at detected the largest shape via many random samples, and the combinatorial complexity of possible samples increases rapidly with the number of required samples.

More specifically, the functions to be implemented are defined in the base class Shape_detection_3::Shape_base:

By default the connected component is detected via the neighbor graph as mentioned above. However, for shapes that admit a faster approach to detect a connected component, the user can provide his/her own implementation to extract the connected component via:

  • Shape_detection_3::Shape_base::connected_component(std::vector<std::size_t>& indices, FT cluster_epsilon): The indices of all supporting points are stored in the vector indices. All points not belonging to the largest cluster of points are removed from the vector indices.

Another optional method can be implemented to provide a helper function providing the shape parameters written to a string:

  • Shape_detection_3::Shape_base::info(): This function returns a string suitable for printing the shape parameters into a log/console. The default solution provides an empty string.

The property maps are used to map the indices to the corresponding points and normals. The following example shows an implementation of a planar shape primitive, which is used by the example Point_set_shape_detection_3/efficient_RANSAC_custom_shape.cpp.


File Point_set_shape_detection_3/efficient_RANSAC_custom_shape.h

#ifndef MY_PLANE_SHAPE_H
#define MY_PLANE_SHAPE_H
#include <CGAL/number_utils.h>
/*
My_Plane derives from Shape_base. The plane is represented by
its normal vector and distance to the origin.
*/
template <class Traits>
class My_Plane : public CGAL::Shape_detection_3::Shape_base<Traits> {
public:
typedef typename Traits::FT FT;
typedef typename Traits::Point_3 Point;
typedef typename Traits::Vector_3 Vector;
public:
My_Plane()
: CGAL::Shape_detection_3::Shape_base<Traits>()
{}
// Computes squared Euclidean distance from query point to the shape.
virtual FT squared_distance(const Point &p) const {
const FT sd = (this->constr_vec(m_point_on_primitive, p)) * m_normal;
return sd * sd;
}
Vector plane_normal() const {
return m_normal;
}
FT d() const {
return m_d;
}
// Returns a string with shape parameters.
virtual std::string info() const {
std::stringstream sstr;
sstr << "Type: plane (" << this->get_x(m_normal) << ", "
<< this->get_y(m_normal) << ", " << this->get_z(m_normal) << ")x - " <<
m_d << " = 0" << " #Pts: " << this->m_indices.size();
return sstr.str();
}
protected:
// Constructs shape based on minimal set of samples from the input data.
virtual void create_shape(const std::vector<std::size_t> &indices) {
const Point p1 = this->point(indices[0]);
const Point p2 = this->point(indices[1]);
const Point p3 = this->point(indices[2]);
m_normal = this->cross_pdct(p1 - p2, p1 - p3);
m_normal = m_normal * (1.0 / sqrt(this->sqlen(m_normal)));
m_d = -(p1[0] * m_normal[0] + p1[1] * m_normal[1] + p1[2] * m_normal[2]);
m_point_on_primitive = p1;
this->m_is_valid = true;
}
// Computes squared Euclidean distance from a set of points.
virtual void squared_distance(const std::vector<std::size_t> &indices,
std::vector<FT> &dists) const {
for (std::size_t i = 0; i < indices.size(); i++) {
const FT sd = (this->point(indices[i])
- m_point_on_primitive) * m_normal;
dists[i] = sd * sd;
}
}
/*
Computes the normal deviation between shape and
a set of points with normals.
*/
virtual void cos_to_normal(const std::vector<std::size_t> &indices,
std::vector<FT> &angles) const {
for (std::size_t i = 0; i < indices.size(); i++)
angles[i] = CGAL::abs(this->normal(indices[i]) * m_normal);
}
// Returns the number of required samples for construction.
virtual std::size_t minimum_sample_size() const {
return 3;
}
private:
Point m_point_on_primitive;
Vector m_normal;
FT m_d;
};
#endif

Plane Regularization

Shape detection applies to man-made scenes or objects such as urban scenes or scans of mechanical parts. Such scenes often contain a wide range of geometric regularities such as parallelism, orthogonality, or symmetry. This package offers a function to reinforce four types of regularities for planar shapes: CGAL::regularize_planes():

  • Planes that are near parallel are made parallel: normal vectors of planes that form angles smaller than a user-defined threshold are made equal.
  • Parallel planes that are near coplanar are made coplanar.
  • Planes that are near orthogonal are made exactly orthogonal.
  • Planes that are near symmetrical with respect to a user-defined axis are made symmetrical.

The user can choose to only regularize one or several of these 4 properties (see reference manual)). The process is greedy and based on a hierarchical decomposition (coplanar clusters are subgroups of parallel clusters which are subgroups of axis-symmetric and orthogonal clusters) as described by Verdie et al. [3]


File Point_set_shape_detection_3/plane_regularization.cpp

#if defined (_MSC_VER) && !defined (_WIN64)
#pragma warning(disable:4244) // boost::number_distance::distance()
// converts 64 to 32 bits integers
#endif
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Point_with_normal_3.h>
#include <CGAL/property_map.h>
#include <iostream>
#include <fstream>
typedef std::pair<Kernel::Point_3, Kernel::Vector_3> Point_with_normal;
typedef std::vector<Point_with_normal> Pwn_vector;
<Kernel, Pwn_vector, Point_map, Normal_map> Traits;
int main(int argc, char** argv)
{
Pwn_vector points;
std::ifstream stream(argc > 1 ? argv[1] : "data/cube.pwn");
if (!stream ||
std::back_inserter(points),
CGAL::parameters::point_map(Point_map()).
normal_map(Normal_map())))
{
std::cerr << "Error: cannot read file cube.pwn" << std::endl;
return EXIT_FAILURE;
}
// Call RANSAC shape detection with planes
Region_growing region_growing;
region_growing.set_input(points);
region_growing.add_shape_factory<Plane>();
region_growing.detect();
Region_growing::Plane_range planes = region_growing.planes();
// Regularize detected planes
Point_map(),
planes,
true, // Regularize parallelism
true, // Regularize orthogonality
false, // Do not regularize coplanarity
true, // Regularize Z-symmetry (default)
10); // 10 degrees of tolerance for parallelism/orthogonality
return EXIT_SUCCESS;
}

Performances

Efficient RANSAC

The running time and detection performance of the efficient RANSAC depend on the chosen parameters. A selective error tolerance parameter leads to higher running times and smaller shapes, as many shape candidates are generated to find the largest shape. We plot the detection performance against the epsilon error tolerance parameter for detecting planes in a complex scene with 5M points, see Figure 72.5. The probability parameter controls the endurance when searching for the largest candidate at each iteration. It barely impacts the number of detected shapes, has a moderate impact on the size of the detected shapes and increases the running times. We plot the performance against the probability parameter, see Figure 72.6.

epsilon_graph.png
Figure 72.5 The graph depicts the number of detected shapes (purple) and the coverage (green), i.e., the ratio assignedPoints / totalPoints, against the epsilon tolerance parameter. A higher value for epsilon, i.e., a more tolerant error, leads to fewer but larger shapes and shorter running times.

prob_graph.png
Figure 72.6 The graph depicts the time, coverage and number of detected primitives against the search endurance parameter, i.e., the probability to miss the largest shape at each iteration. The number of shapes is stable and the coverage increases when the probability is lowered. The running times increase significantly as many more candidates are generated during each iteration of the algorithm.

Region Growing

Region growing iterates through every single point in the input, which makes it usually slower (although more robust). The main parameter that has an effect on running times is the epsilon used for clustering: internally, range queries are perform using spheres with radius cluster_epsilon. Using larger values means using larger spheres which are more computationally demanding.

We plot again the detection performance against the epsilon parameter for detecting planes, this time in a scene of 500k points, see figure Figure 72.7.

epsilon_graph_2.png
Figure 72.7 The graph depicts the number of detected shapes (purple) and the coverage (green), i.e., the ratio assignedPoints / totalPoints, against the epsilon parameter (for simplicity, we use the same value for epsilon and cluster_epsilon). A higher value for cluster_epsilon, i.e., larger neighborhood spheres, leads to longer computation and to larger shapes.

History

The efficient RANSAC implementation was developed by Sven Oesau based on a prototype version by Yannick Verdie, with the help of Clément Jamin and under the supervision of Pierre Alliez. Plane regularization and region growing were added by Simon Giraudot based on prototype versions developed by Florent Lafarge.