\( \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 5.0.2 - Shape Detection
User Manual

Authors
Sven Oesau, Yannick Verdie, Clément Jamin, Pierre Alliez, Florent Lafarge, Simon Giraudot, Thien Hoang, and Dmitry Anisimov

Introduction

This CGAL component implements two algorithms for shape detection:

  • the Efficient RANSAC (RANdom SAmple Consensus) method, contributed by Schnabel et al. [2];
  • the Region Growing method, contributed by Lafarge and Mallet [1].

Efficient RANSAC

From an unstructured point set with unoriented normals, this algorithm detects a set of shapes (see Figure Figure 77.1). Five types of primitive shapes are provided by this package: plane, sphere, cylinder, cone, and torus. Other primitive shapes can be easily added by the user (see Section Custom Shapes).

overview2.png
Figure 77.1 Input and output of the Efficient RANSAC method. (a) Input point set. (b) Point set depicted with one color per detected shape.

This method takes as input a point set with unoriented normals and provides as output a set of detected shapes with associated input points. The output of the algorithm 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.

The shapes are detected via a RANSAC-type approach, that is 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 the highest number of inliers, referred to as the 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 minimum 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 is 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 threshold. 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 the 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).

Parameters

The algorithm has five parameters:

  • 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\). The parameter epsilon defines the absolute maximum tolerance Euclidean distance between a point and a shape. A high value of epsilon leads to the detection of fewer large shapes and hence a less detailed detection. A low value of 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 enables to detect more primitive shapes that approximate some of the free-form surfaces. The impact of this parameter is depicted by Figure Figure 77.2. Its impact on performance is evaluated in Section Performance.

epsilon_variation2.png
Figure 77.2 Impact of the epsilon parameter over the levels of detail 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 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. The 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 the parameter cluster_epsilon is depicted in Figure Figure 77.3.

varying_connectivity.png
Figure 77.3 The 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 of cluster_epsilon leads to detecting a single planar shape. (b) A moderate value of cluster_epsilon yields the detection of four squares. Notice that a few points within the squares are not detected as not connected. (c) A small value of cluster_epsilon leads to over-segmentation.

  • min_points: The 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 this parameter is not strict: depending on the chosen probability, shapes may be extracted with a number of points lower than the specified 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 time due to a higher search endurance.

Examples

The main class Shape_detection::Efficient_RANSAC takes a template parameter Shape_detection::Efficient_RANSAC_traits that defines the geometric types and input format. Property maps provide a means to interface with the user-specific data structures. The first parameter of the Shape_detection::Efficient_RANSAC_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::Efficient_RANSAC_traits class. The concept behind property maps is detailed in Manual CGAL and Property 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 Plane Detection

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


File Shape_detection/efficient_RANSAC_basic.cpp

#if defined (_MSC_VER) && !defined (_WIN64)
#pragma warning(disable:4244) // boost::number_distance::distance()
// converts 64 to 32 bits integers
#endif
#include <fstream>
#include <iostream>
#include <CGAL/property_map.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Point_with_normal_3.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Shape_detection/Efficient_RANSAC.h>
// Type declarations.
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) {
std::cout << "Efficient RANSAC" << std::endl;
const char* filename = (argc > 1) ? argv[1] : "data/cube.pwn";
// Points with normals.
Pwn_vector points;
// Load point set from a file.
std::ifstream stream(filename);
if (!stream ||
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;
}
// Instantiate shape detection engine.
Efficient_ransac ransac;
// Provide input data.
ransac.set_input(points);
// Register planar shapes via template method.
ransac.add_shape_factory<Plane>();
// Detect registered shapes with default parameters.
ransac.detect();
// Print number of detected shapes.
std::cout << ransac.shapes().end() - ransac.shapes().begin()
<< " shapes detected." << std::endl;
return EXIT_SUCCESS;
}

Plane Detection With Callback

The Efficient RANSAC class provides a callback mechanism that enables the user to track the progress of the algorithm. 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 Shape_detection/efficient_RANSAC_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 <fstream>
#include <iostream>
#include <CGAL/Timer.h>
#include <CGAL/property_map.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Point_with_normal_3.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Shape_detection/Efficient_RANSAC.h>
// Type declarations.
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;
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 reached, interrupt the algorithm.
if (timer.time() > limit) {
std::cerr << "Algorithm takes too long, exiting ("
<< 100.0 * advancement << "% done)" << std::endl;
return false;
}
return true;
}
};
int main (int argc, char** argv) {
std::cout << "Efficient RANSAC" << std::endl;
const char* filename = (argc > 1) ? argv[1] : "data/cube.pwn";
Pwn_vector points;
std::ifstream stream(filename);
if (!stream ||
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;
}
Efficient_ransac ransac;
ransac.set_input(points);
ransac.add_shape_factory<Plane>();
// Create callback that interrupts the algorithm
// if it takes more than half a second.
Timeout_callback timeout_callback(0.5);
// Detect registered shapes with the default parameters.
ransac.detect(Efficient_ransac::Parameters(), timeout_callback);
return EXIT_SUCCESS;
}

Setting Parameters And Using Different Shape Types

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

Basic information of the detected shapes is written to the standard output: if the shape 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 shapes.


File Shape_detection/efficient_RANSAC_with_parameters.cpp

#include <fstream>
#include <iostream>
#include <CGAL/property_map.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Point_with_normal_3.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Shape_detection/Efficient_RANSAC.h>
// 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;
<Kernel, Pwn_vector, Point_map, Normal_map> Traits;
int main(int argc, char** argv) {
// Points with normals.
Pwn_vector points;
// Load point set from a file.
std::ifstream stream((argc > 1) ? argv[1] : "data/cube.pwn");
if (!stream ||
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;
// Instantiate shape detection engine.
Efficient_ransac ransac;
// Provide 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>();
// Set parameters for shape detection.
Efficient_ransac::Parameters parameters;
// Set probability to miss the largest primitive at each iteration.
parameters.probability = 0.05;
// Detect shapes with at least 500 points.
parameters.min_points = 200;
// Set maximum Euclidean distance between a point and a shape.
parameters.epsilon = 0.002;
// Set maximum Euclidean distance between points to be clustered.
parameters.cluster_epsilon = 0.01;
// Set maximum normal deviation.
// 0.9 < dot(surface_normal, point_normal);
parameters.normal_threshold = 0.9;
// Detect shapes.
ransac.detect(parameters);
// Print 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 the 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 the 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 {
// Print the parameters of the detected shape.
// This function is available for any type of shape.
std::cout << (*it)->info() << std::endl;
}
// Proceed with the 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 Shape_detection/efficient_RANSAC_with_point_access.cpp

#include <fstream>
#include <iostream>
#include <CGAL/Timer.h>
#include <CGAL/number_utils.h>
#include <CGAL/property_map.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Shape_detection/Efficient_RANSAC.h>
// 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;
<Kernel, Pwn_vector, Point_map, Normal_map> Traits;
int main(int argc, char** argv) {
// Points with normals.
Pwn_vector points;
// Load point set from a file.
std::ifstream stream((argc > 1) ? argv[1] : "data/cube.pwn");
if (!stream ||
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;
}
// Instantiate shape detection engine.
Efficient_ransac ransac;
// Provide input data.
ransac.set_input(points);
// Register detection of planes.
ransac.add_shape_factory<Plane>();
// Measure time before setting up the shape detection.
CGAL::Timer time;
time.start();
// Build internal data structures.
ransac.preprocess();
// Measure time after preprocessing.
time.stop();
std::cout << "preprocessing took: " << time.time() * 1000 << "ms" << std::endl;
// Perform detection several times and choose result with the highest coverage.
Efficient_ransac::Shape_range shapes = ransac.shapes();
FT best_coverage = 0;
for (std::size_t i = 0; i < 3; ++i) {
// Reset timer.
time.reset();
time.start();
// Detect shapes.
ransac.detect();
// Measure 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());
// Print number of assigned shapes and unassigned 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 the 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;
// Use Shape_base::info() to print the parameters of the detected shape.
std::cout << (*it)->info();
// Sums distances of points to the detected shapes.
FT sum_distances = 0;
// Iterate 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()) {
// Retrieve 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));
// Proceed with the next point.
index_it++;
}
// Compute and print the average distance.
FT average_distance = sum_distances / shape->indices_of_assigned_points().size();
std::cout << " average distance: " << average_distance << std::endl;
// Proceed with the next detected shape.
it++;
}
return EXIT_SUCCESS;
}

Custom Shapes

Other shape types can be detected by implementing a shape class derived from the class Shape_detection::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 to the shape, and compute the normal deviation between a query point with the 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 detecting 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::Shape_base:

The last two functions are used to determine the number of inlier points to the shape. They compute respectively the squared distance from a set of points to the shape, and the dot product between the point normals and the normals at the shape for the closest points on the shape.

The access to the actual point and normal data is carried out via Shape_detection::Shape_base::point(std::size_t index) and Shape_detection::Shape_base::normal(std::size_t index) (see the example below). The resulting squared distance/dot product is stored in the vector provided as the first argument.

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::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 that do not belong 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::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 header shows an implementation of a planar shape primitive, which is used by the example Shape_detection/efficient_RANSAC_with_custom_shape.cpp.


File Shape_detection/include/efficient_RANSAC_with_custom_shape.h

#ifndef MY_PLANE_SHAPE_H
#define MY_PLANE_SHAPE_H
#include <CGAL/number_utils.h>
#include <CGAL/Shape_detection/Efficient_RANSAC.h>
// My_Plane is derived 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::Shape_base<Traits> {
public:
typedef typename Traits::FT FT; // number type
typedef typename Traits::Point_3 Point; // point type
typedef typename Traits::Vector_3 Vector; // vector type
My_Plane() :
CGAL::Shape_detection::Shape_base<Traits>()
{ }
// Compute 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;
}
// Return 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:
// Construct shape base on a 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;
}
// Compute 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;
}
}
// Compute the normal deviation between a 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);
}
// Return 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 // MY_PLANE_SHAPE_H

Performance

The running time and detection performance of the Efficient RANSAC depend on the chosen parameters. A selective error tolerance parameter leads to higher running time and fewer 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 Figure 77.4). 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 time. We plot the performance against the probability parameter (see Figure Figure 77.5).

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

prob_graph.png
Figure 77.5 The graph depicts the time, coverage, and the number of detected primitives against the search endurance parameter, that is 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 time increases significantly as many more candidates are generated during each iteration of the algorithm.

Region Growing

This shape detection component is based on the region growing algorithm applied to a set of user-specified items. Shapes are detected by growing regions from seed items, where each region is created as follows:

  1. Pick the next available seed item;
  2. Find its neighbors in the data set;
  3. Include those neighbors, which satisfy the region requirements;
  4. Repeat the procedure for all included neighbors;
  5. If no further neighbor satisfies the requirements, start a new region.

Together with the generic algorithm's implementation CGAL::Shape_detection::Region_growing, three particular instances of this algorithm are provided:

Other instances can be easily added by the user, as explained below.

Framework

The main class CGAL::Shape_detection::Region_growing is parameterized by

  • InputRange that stores a range of user-defined input items;
  • NeighborQuery that provides the means for accessing neighbors of an item;
  • RegionType that provides the means for validating regions;
  • SeedMap that defines the seeding order of items.

Using this generic framework, users can grow any type of regions on a set of arbitrary items with their own propagation and seeding conditions (see an example).

Neighborhood

The concept NeighborQuery provides the means for accessing neighbors of an item. To create a model that respects this concept, the user has to provide an overload of the operator:

Regions

The concept RegionType provides the means for validating regions. In fact, a model of this concept maintains a description of the region type that is used in region growing. To create a model that respect this concept, three functions have to be defined:

  • RegionType::is_part_of_region() This function checks if an item satisfies all necessary region requirements and can be added to a region. It is called per item.
  • RegionType::is_valid_region() This function checks if a region satisfies all necessary region requirements. It is called per region.
  • RegionType::update() This utility function enables to update any information, which is maintained with the region.

Seeding

The SeedMap property map, provided as an optional parameter to the main class, enables to define the seeding order of items that is which items are used first to grow regions from. Such items are referred to as seed items. The SeedMap maps the index of an item to its order number in the overall region growing processing queue. The default map is the identity one that is the seed index of the item equals to the item's index in the input_range. If it maps to std::size_t(-1), then the corresponding item is skipped.

Examples

This toy example shows how to define one's own NeighborQuery and RegionType classes, which are used to parameterize the CGAL::Shape_detection::Region_growing. It also shows how to skip unnecessary items and change their default seeding order.

We choose a simple custom object item. We define four such objects, where for each object, we manually store indices of its neighbors. The operator NeighborQuery::operator()() does nothing but accessing these neighbors. The RegionType class defines the three necessary functions:

We also define a SeedMap, such that the second object is handled first, while the first object follows. Moreover, the last object is always skipped. Notice that in this example, the container with objects is std::list, which is not a random access. In general, it is much slower, but illustrates that the generic region growing algorithm can be applied to any type of a container.

The result of using these classes with the region growing main class is that the first two objects form the first region, the third object forms the second region, and the last object is skipped.


File Shape_detection/region_growing_with_custom_classes.cpp

// STL includes.
#include <map>
#include <vector>
#include <string>
#include <iostream>
#include <iterator>
// CGAL includes.
#include <CGAL/Shape_detection/Region_growing/Region_growing.h>
// Custom Neighbor_query, Region_type, and Seed_map classes for region growing.
namespace Custom {
// An object that stores indices of all its neighbors.
struct Object {
std::vector<std::size_t> neighbors;
};
// A range of objects.
using Objects = std::vector<Object>;
// The Neighbor_query functor that accesses neighbors stored in
// the object struct above.
class Neighbor_query {
const Objects& m_objects;
public:
Neighbor_query(const Objects& objects) :
m_objects(objects)
{ }
void operator()(
const std::size_t query_index,
std::vector<std::size_t>& neighbors) const {
std::size_t i = 0;
for (const auto& object : m_objects) {
if (i == query_index) {
neighbors = object.neighbors;
return;
} ++i;
}
}
};
// The Region_type class, where the function is_part_of_region() verifies
// a very specific condition that the first and second objects in the
// range are in fact neighbors; is_valid_region() function always
// returns true after the first call to the update() function.
// These are the only functions that have to be defined.
class Region_type {
bool m_is_valid = false;
public:
Region_type() { }
bool is_part_of_region(
const std::size_t,
const std::size_t query_index,
const std::vector<std::size_t>& region) const {
if (region.size() == 0)
return false;
const std::size_t index = region[0];
if (index == 0 && query_index == 1) return true;
if (query_index == 0 && index == 1) return true;
return false;
}
inline bool is_valid_region(const std::vector<std::size_t>&) const {
return m_is_valid;
}
void update(const std::vector<std::size_t>&) {
m_is_valid = true;
}
};
// The SeedMap class that uses the m_objects_map to define
// the seeding order of objects.
class Seed_map {
public:
using key_type = std::size_t;
using value_type = std::size_t;
using category = boost::lvalue_property_map_tag;
Seed_map(const std::map<std::size_t, std::size_t>& objects_map) :
m_objects_map(objects_map)
{ }
value_type operator[](const key_type key) const {
return m_objects_map.find(key)->second;
}
friend value_type get(
const Seed_map& seed_map,
const key_type key) {
return seed_map[key];
}
private:
const std::map<std::size_t, std::size_t>& m_objects_map;
};
} // namespace Custom
// Type declarations.
using Object = Custom::Object;
using Objects = Custom::Objects;
using Neighbor_query = Custom::Neighbor_query;
using Region_type = Custom::Region_type;
using Seed_map = Custom::Seed_map;
using Region = std::vector<std::size_t>;
using Regions = std::vector<Region>;
using Region_growing =
int main() {
std::cout << std::endl <<
"region_growing_with_custom_classes example started"
<< std::endl << std::endl;
// Define a range of objects, where the first two objects form
// the first region, while the third object forms the second region.
// Note that Objects is a random access container here, however the
// same algorithm/example can work with other containers, e.g. std::list.
Objects objects;
// Region 1.
Object object1;
object1.neighbors.resize(1, 1);
objects.push_back(object1);
Object object2;
object2.neighbors.resize(1, 0);
objects.push_back(object2);
// Region 2.
Object object3;
objects.push_back(object3);
// Extra object to skip.
Object object4;
objects.push_back(object4);
// Create instances of the classes Neighbor_query and Region_type.
Neighbor_query neighbor_query = Neighbor_query(objects);
Region_type region_type = Region_type();
// Create a seed map.
std::map<std::size_t, std::size_t> objects_map;
objects_map[0] = 1; // the order is swapped with the next object
objects_map[1] = 0;
objects_map[2] = 2; // the default order
objects_map[3] = std::size_t(-1); // skip this object
const Seed_map seed_map(objects_map);
// Create an instance of the region growing class.
Region_growing region_growing(
objects, neighbor_query, region_type, seed_map);
// Run the algorithm.
Regions regions;
region_growing.detect(std::back_inserter(regions));
// Print the number of found regions. It must be two regions.
std::cout << "* " << regions.size() <<
" regions have been found among " << objects.size() << " objects"
<< std::endl;
std::cout << std::endl <<
"region_growing_with_custom_classes example finished"
<< std::endl << std::endl;
return EXIT_SUCCESS;
}

Point Set

If one wants to detect lines (see 2D Example)

region_growing_on_point_set_2.png
Figure 77.6 A 2D point set depicted with one color per detected line.

or planes (see 3D Example)

region_growing_on_point_set_3.png
Figure 77.7 A 3D point set depicted with one color per detected plane.

in a 2D or 3D point set respectively, this CGAL component provides the corresponding models of the concepts NeighborQuery and RegionType. In particular, it provides two different ways to define neighbors of a point:

The component also provides

The program associates all points from a region to the best-fit hyperplane (2D line or 3D plane) and controls the quality of this fit.

The quality of region growing in a point set (2D or 3D) can be improved by slightly sacrificing the running time. To achieve this, one can sort indices of input points with respect to some quality criteria. These quality criteria can be included through the SeedMap input parameter. We provide a quality sorting both for 2D and 3D points:

Parameters

The classes in Section Region Growing On Point Set depend on a few parameters that should be defined by the user. They also have default values, but these values do not necessarily guarantee to produce pleasant results.

The NeighborQuery related classes depend on the following parameters:

The right choice of sphere_radius or k parameters plays an important role in producing a good result. For example, if we consider the fuzzy sphere neighborhood, when sphere_radius is too large, we have fewer regions, and the details are not clearly separated. Meanwhile, if sphere_radius is too small, we produce more regions, and the point set may be over-segmented. Consider a 2D map of an intersection of streets in a city as in Figure Figure 77.8. Each region is painted with a unique color. As sphere_radius increases, the details become less clear. When sphere_radius = 0.3 (c), the best visual result is produced.

sphere_radius_parameter_2D.png
Figure 77.8 (a) Input 2D point set; (b) 17 regions are found when sphere_radius = 0.1; (c) 8 regions are found when sphere_radius = 0.3; (d) 4 regions are found when sphere_radius = 1.2.

The RegionType related classes depend on the following parameters:

  • distance_threshold - the maximum distance from a point to a line/plane;
  • angle_threshold - the maximum accepted angle between the normal associated with a point and the normal of a line/plane;
  • min_region_size - the minimum number of points a region must have.

The first two parameters are used by the functions RegionType::is_part_of_region() and RegionType::update(), while the third parameter is used by the function RegionType::is_valid_region() explained in Section Framework Region Type.

The right choice of distance_threshold and angle_threshold parameters is also very important. For example, Figure Figure 77.9 shows that the roof top of the house can be distinguished as two planes (painted in blue and dark red) when angle_threshold is strict enough (c), or it can be recognized as only one plane (painted in pale yellow) in the other case (b).

angle_threshold_parameter_3D.png
Figure 77.9 (a) Input 3D point cloud; (b) Result when angle_threshold = 60 degrees; (c) Result when angle_threshold = 25 degrees.

Examples

Typical usage of the Region Growing component consists of five steps:

  1. Define an input range with points;
  2. Create instances of the classes NeighborQuery and RegionType with the proper parameters;
  3. Create an instance of the class CGAL::Shape_detection::Region_growing;
  4. Detect;
  5. Postprocess.

Given a 2D point set, we detect 2D lines using the fuzzy sphere neighborhood. We then color all points from the found regions and save them in a file (see Figure Figure 77.6).

The points with assigned to them normal vectors are stored in std::vector and the used Kernel is CGAL::Simple_cartesian, where the number type is double.


File Shape_detection/region_growing_on_point_set_2.cpp

// STL includes.
#include <string>
#include <vector>
#include <utility>
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <iterator>
// CGAL includes.
#include <CGAL/array.h>
#include <CGAL/property_map.h>
#include <CGAL/IO/write_ply_points.h>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Shape_detection/Region_growing/Region_growing.h>
#include <CGAL/Shape_detection/Region_growing/Region_growing_on_point_set.h>
// Type declarations.
using FT = typename Kernel::FT;
using Point_2 = typename Kernel::Point_2;
using Point_3 = typename Kernel::Point_3;
using Vector_2 = typename Kernel::Vector_2;
using Point_with_normal = std::pair<Point_2, Vector_2>;
using Input_range = std::vector<Point_with_normal>;
using Region = std::vector<std::size_t>;
using Regions = std::vector<Region>;
using Color = std::array<unsigned char, 3>;
using Point_with_color = std::pair<Point_3, Color>;
using Pwc_vector = std::vector<Point_with_color>;
// Define how a color should be stored.
namespace CGAL {
template<class F>
struct Output_rep< ::Color, F > {
const ::Color& c;
static const bool is_specialized = true;
Output_rep(const ::Color& c) : c(c) { }
std::ostream& operator()(std::ostream& out) const {
if (is_ascii(out))
out << int(c[0]) << " " << int(c[1]) << " " << int(c[2]);
else
out.write(reinterpret_cast<const char*>(&c), sizeof(c));
return out;
}
};
} // namespace CGAL
int main(int argc, char *argv[]) {
std::cout << std::endl <<
"region_growing_on_point_set_2 example started"
<< std::endl << std::endl;
std::cout <<
"Note: if 0 points are loaded, please specify the path to the file data/point_set_2.xyz by hand!"
<< std::endl << std::endl;
// Load xyz data either from a local folder or a user-provided file.
std::ifstream in(argc > 1 ? argv[1] : "data/point_set_2.xyz");
Input_range input_range;
FT a, b, c, d, e, f;
while (in >> a >> b >> c >> d >> e >> f)
input_range.push_back(
std::make_pair(Point_2(a, b), Vector_2(d, e)));
in.close();
std::cout <<
"* loaded "
<< input_range.size() <<
" points with normals"
<< std::endl;
// Default parameter values for the data file point_set_2.xyz.
const FT search_sphere_radius = FT(5);
const FT max_distance_to_line = FT(45) / FT(10);
const FT max_accepted_angle = FT(45);
const std::size_t min_region_size = 5;
// Create instances of the classes Neighbor_query and Region_type.
Neighbor_query neighbor_query(
input_range,
search_sphere_radius);
Region_type region_type(
input_range,
max_distance_to_line, max_accepted_angle, min_region_size);
// Create an instance of the region growing class.
Region_growing region_growing(
input_range, neighbor_query, region_type);
// Run the algorithm.
Regions regions;
region_growing.detect(std::back_inserter(regions));
// Print the number of found regions.
std::cout << "* " << regions.size() <<
" regions have been found"
<< std::endl;
Pwc_vector pwc;
srand(static_cast<unsigned int>(time(NULL)));
// Iterate through all regions.
for (const auto& region : regions) {
// Generate a random color.
const Color color =
static_cast<unsigned char>(rand() % 256),
static_cast<unsigned char>(rand() % 256),
static_cast<unsigned char>(rand() % 256));
// Iterate through all region items.
for (const auto index : region) {
const auto& key = *(input_range.begin() + index);
const Point_2& point = get(Point_map(), key);
pwc.push_back(std::make_pair(Point_3(point.x(), point.y(), 0), color));
}
}
// Save the result to a file in the user-provided path if any.
if (argc > 2) {
const std::string path = argv[2];
const std::string fullpath = path + "regions_point_set_2.ply";
std::ofstream out(fullpath);
out, pwc,
CGAL::make_ply_point_writer(PLY_Point_map()),
std::make_tuple(
PLY_Color_map(),
std::cout <<
"* found regions are saved in "
<< fullpath << std::endl;
out.close();
}
std::cout << std::endl <<
"region_growing_on_point_set_2 example finished"
<< std::endl << std::endl;
return EXIT_SUCCESS;
}

If we are given a 3D point set, then the example below shows how to detect 3D planes using the K nearest neighbors search. We color all points from the found regions and save them in a file (see Figure Figure 77.7). The example also shows how to retrieve all points, which are not assigned to any region, and how to use a custom output iterator.

The point set with associated normals is stored in CGAL::Point_set_3 and the used Kernel is CGAL::Exact_predicates_inexact_constructions_kernel.


File Shape_detection/region_growing_on_point_set_3.cpp

// STL includes.
#include <string>
#include <vector>
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <iterator>
// Boost includes.
#include <boost/function_output_iterator.hpp>
// CGAL includes.
#include <CGAL/Timer.h>
#include <CGAL/Random.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Point_set_3.h>
#include <CGAL/Point_set_3/IO.h>
#include <CGAL/Shape_detection/Region_growing/Region_growing.h>
#include <CGAL/Shape_detection/Region_growing/Region_growing_on_point_set.h>
// Type declarations.
using FT = typename Kernel::FT;
using Point_3 = typename Kernel::Point_3;
using Vector_3 = typename Kernel::Vector_3;
using Input_range = CGAL::Point_set_3<Point_3>;
using Point_map = typename Input_range::Point_map;
using Normal_map = typename Input_range::Vector_map;
using Indices = std::vector<std::size_t>;
using Output_range = CGAL::Point_set_3<Point_3>;
using Points_3 = std::vector<Point_3>;
// Define an insert iterator.
struct Insert_point_colored_by_region_index {
using argument_type = Indices;
using result_type = void;
using Color_map =
typename Output_range:: template Property_map<unsigned char>;
const Input_range& m_input_range;
const Point_map m_point_map;
Output_range& m_output_range;
std::size_t& m_number_of_regions;
Color_map m_red, m_green, m_blue;
Insert_point_colored_by_region_index(
const Input_range& input_range,
const Point_map point_map,
Output_range& output_range,
std::size_t& number_of_regions) :
m_input_range(input_range),
m_point_map(point_map),
m_output_range(output_range),
m_number_of_regions(number_of_regions) {
m_red =
m_output_range.template add_property_map<unsigned char>("red", 0).first;
m_green =
m_output_range.template add_property_map<unsigned char>("green", 0).first;
m_blue =
m_output_range.template add_property_map<unsigned char>("blue", 0).first;
}
result_type operator()(const argument_type& region) {
CGAL::Random rand(static_cast<unsigned int>(m_number_of_regions));
const unsigned char r =
static_cast<unsigned char>(64 + rand.get_int(0, 192));
const unsigned char g =
static_cast<unsigned char>(64 + rand.get_int(0, 192));
const unsigned char b =
static_cast<unsigned char>(64 + rand.get_int(0, 192));
for (const std::size_t index : region) {
const auto& key = *(m_input_range.begin() + index);
const Point_3& point = get(m_point_map, key);
const auto it = m_output_range.insert(point);
m_red[*it] = r;
m_green[*it] = g;
m_blue[*it] = b;
}
++m_number_of_regions;
}
}; // Insert_point_colored_by_region_index
int main(int argc, char *argv[]) {
std::cout << std::endl <<
"region_growing_on_point_set_3 example started"
<< std::endl << std::endl;
std::cout <<
"Note: if 0 points are loaded, please specify the path to the file data/point_set_3.xyz by hand!"
<< std::endl << std::endl;
// Load xyz data either from a local folder or a user-provided file.
std::ifstream in(argc > 1 ? argv[1] : "data/point_set_3.xyz");
const bool with_normal_map = true;
Input_range input_range(with_normal_map);
in >> input_range;
in.close();
std::cout <<
"* loaded "
<< input_range.size() <<
" points with normals"
<< std::endl;
// Default parameter values for the data file point_set_3.xyz.
const std::size_t k = 12;
const FT max_distance_to_plane = FT(2);
const FT max_accepted_angle = FT(20);
const std::size_t min_region_size = 50;
// Create instances of the classes Neighbor_query and Region_type.
Neighbor_query neighbor_query(
input_range,
k,
input_range.point_map());
Region_type region_type(
input_range,
max_distance_to_plane, max_accepted_angle, min_region_size,
input_range.point_map(), input_range.normal_map());
// Create an instance of the region growing class.
Region_growing region_growing(
input_range, neighbor_query, region_type);
// Run the algorithm.
Output_range output_range;
std::size_t number_of_regions = 0;
Insert_point_colored_by_region_index inserter(
input_range, input_range.point_map(),
output_range, number_of_regions);
CGAL::Timer timer;
timer.start();
region_growing.detect(
boost::make_function_output_iterator(inserter));
timer.stop();
// Print the number of found regions.
std::cout << "* " << number_of_regions <<
" regions have been found in " << timer.time() << " seconds"
<< std::endl;
// Save the result to a file in the user-provided path if any.
if (argc > 2) {
const std::string path = argv[2];
const std::string fullpath = path + "regions_point_set_3.ply";
std::ofstream out(fullpath);
out << output_range;
std::cout << "* found regions are saved in " << fullpath << std::endl;
out.close();
}
// Get all unassigned items.
Indices unassigned_items;
region_growing.unassigned_items(std::back_inserter(unassigned_items));
// Print the number of unassigned items.
std::cout << "* " << unassigned_items.size() <<
" points do not belong to any region"
<< std::endl;
// Store all unassigned points.
Points_3 unassigned_points;
unassigned_points.reserve(unassigned_items.size());
for (const auto index : unassigned_items) {
const auto& key = *(input_range.begin() + index);
const Point_3& point = get(input_range.point_map(), key);
unassigned_points.push_back(point);
}
std::cout << "* " << unassigned_points.size() <<
" unassigned points are stored"
<< std::endl;
std::cout << std::endl <<
"region_growing_on_point_set_3 example finished"
<< std::endl << std::endl;
return EXIT_SUCCESS;
}

Performance

The main parameter that affects the region growing algorithm on a point set is the neighborhood size at each retrieval (sphere_radius or k). Larger neighbor lists are often followed by a smaller number of regions, larger coverage (the ratio between the number of points assigned to regions and the total number of input points), and longer running time. For example, for a test of about 70k 2D points with the fuzzy sphere neighborhood, the following table is produced:

sphere_radius Time (in seconds) Number of regions Number of assigned points
1 0.138831 794 4483
3 0.069098 3063 63038
6 0.077703 2508 64906
9 0.093415 2302 65334

If the neighborhood size is set too low, some points might be isolated, the region size would not reach a critical mass and so will be discarded. This does not only cause the latency in the program, but also reduces the coverage value, as can be seen when the sphere_radius = 1. A typical time measure for a 3D point set with the K nearest neighborhood and well-defined parameters can be found in the following table:

Number of points Time (in seconds)
300k 0.761617
600k 1.68735
900k 2.80346
1200k 4.06246

Polygon Mesh

If one wants to detect planes on a polygon mesh, this CGAL component provides the corresponding models of the concepts NeighborQuery and RegionType. In particular, it has

This component accepts any model of the concept FaceListGraph as a polygon mesh. A picture below gives an example of the region growing algorithm for detecting 3D planes on CGAL::Surface_mesh.

region_growing_on_polygon_mesh.png
Figure 77.10 A surface mesh depicted with one color per detected plane.

The quality of region growing on a polygon mesh can be improved by slightly sacrificing the running time. To achieve this, one can sort indices of input faces with respect to some quality criteria. These quality criteria can be included in region growing through the SeedMap input parameter. We provide such a quality sorting:

Parameters

The NeighborQuery related class does not require any parameters, because edge-adjacent faces are found using the internal face graph connectivity, while the RegionType related class depends on three parameters:

  • distance_threshold - the maximum distance from the furthest face vertex to a plane;
  • angle_threshold - the maximum accepted angle between the face normal and the normal of a plane;
  • min_region_size - the minimum number of mesh faces a region must have.

The first two parameters are used by the functions RegionType::is_part_of_region() and RegionType::update(), while the third parameter is used by the function RegionType::is_valid_region() explained in Section Framework Regions. The right choice of these parameters is as important as the one explained in Section Parameters For Region Growing On Point Set.

Examples

In the example below, we show how to use region growing to detect planes on a polygon mesh that can be either stored as CGAL::Surface_mesh or CGAL::Polyhedron_3. If it is a surface mesh, this example also provides a way to save the result in a file (see Figure Figure 77.10). The used Kernel here is CGAL::Exact_predicates_exact_constructions_kernel. Though this exact kernel provides better quality results, please note that it may significantly slow down the execution of the program, so if you need a faster version of region growing, use a floating type based kernel.

We can improve the quality of region growing by providing a different seeding order (analogously to Point Sets) that is why in this example we also sort indices of input faces using the CGAL::Shape_detection::Polygon_mesh::Least_squares_plane_fit_sorting and only then detect regions.


File Shape_detection/region_growing_on_polygon_mesh.cpp

// STL includes.
#include <string>
#include <vector>
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <iterator>
// CGAL includes.
#include <CGAL/memory.h>
#include <CGAL/IO/Color.h>
#include <CGAL/Iterator_range.h>
#include <CGAL/HalfedgeDS_vector.h>
#include <CGAL/Exact_predicates_exact_constructions_kernel.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/Polyhedron_3.h>
#include <CGAL/Shape_detection/Region_growing/Region_growing.h>
#include <CGAL/Shape_detection/Region_growing/Region_growing_on_polygon_mesh.h>
// Type declarations.
using FT = typename Kernel::FT;
using Point_3 = typename Kernel::Point_3;
using Color = CGAL::Color;
// Choose the type of a container for a polygon mesh.
#define USE_SURFACE_MESH
#if defined(USE_SURFACE_MESH)
using Polygon_mesh = CGAL::Surface_mesh<Point_3>;
using Face_range = typename Polygon_mesh::Face_range;
#else
#endif
using Region = std::vector<std::size_t>;
using Regions = std::vector<Region>;
using Vertex_to_point_map = typename Region_type::Vertex_to_point_map;
int main(int argc, char *argv[]) {
std::cout << std::endl <<
"region_growing_on_polygon_mesh example started"
<< std::endl << std::endl;
// Load off data either from a local folder or a user-provided file.
std::ifstream in(argc > 1 ? argv[1] : "data/polygon_mesh.off");
Polygon_mesh polygon_mesh;
in >> polygon_mesh;
in.close();
const Face_range face_range = faces(polygon_mesh);
std::cout <<
"* polygon mesh with "
<< face_range.size() <<
" faces is loaded"
<< std::endl;
// Default parameter values for the data file polygon_mesh.off.
const FT max_distance_to_plane = FT(1);
const FT max_accepted_angle = FT(45);
const std::size_t min_region_size = 5;
// Create instances of the classes Neighbor_query and Region_type.
Neighbor_query neighbor_query(polygon_mesh);
const Vertex_to_point_map vertex_to_point_map(
get(CGAL::vertex_point, polygon_mesh));
Region_type region_type(
polygon_mesh,
max_distance_to_plane, max_accepted_angle, min_region_size,
vertex_to_point_map);
// Sort face indices.
Sorting sorting(
polygon_mesh, neighbor_query,
vertex_to_point_map);
sorting.sort();
// Create an instance of the region growing class.
Region_growing region_growing(
face_range, neighbor_query, region_type,
sorting.seed_map());
// Run the algorithm.
Regions regions;
region_growing.detect(std::back_inserter(regions));
// Print the number of found regions.
std::cout << "* " << regions.size() <<
" regions have been found"
<< std::endl;
// Save the result in a file only if it is stored in CGAL::Surface_mesh.
#if defined(USE_SURFACE_MESH)
using Face_index = typename Polygon_mesh::Face_index;
// Save the result to a file in the user-provided path if any.
srand(static_cast<unsigned int>(time(NULL)));
if (argc > 2) {
bool created;
typename Polygon_mesh::template Property_map<Face_index, Color> face_color;
boost::tie(face_color, created) =
polygon_mesh.template add_property_map<Face_index, Color>(
"f:color", Color(0, 0, 0));
if (!created) {
std::cout << std::endl <<
"region_growing_on_polygon_mesh example finished"
<< std::endl << std::endl;
return EXIT_FAILURE;
}
const std::string path = argv[2];
const std::string fullpath = path + "regions_polygon_mesh.off";
std::ofstream out(fullpath);
// Iterate through all regions.
for (const auto& region : regions) {
// Generate a random color.
const Color color(
static_cast<unsigned char>(rand() % 256),
static_cast<unsigned char>(rand() % 256),
static_cast<unsigned char>(rand() % 256));
// Iterate through all region items.
using size_type = typename Polygon_mesh::size_type;
for (const auto index : region)
face_color[Face_index(static_cast<size_type>(index))] = color;
}
out << polygon_mesh;
out.close();
std::cout <<
"* polygon mesh is saved in "
<< fullpath << std::endl;
}
#endif
std::cout << std::endl <<
"region_growing_on_polygon_mesh example finished"
<< std::endl << std::endl;
return EXIT_SUCCESS;
}

Performance

Since accessing neighbors of a face in a polygon mesh is fast, performance of the region growing on a polygon mesh mostly depends on the RegionType model's implementation, which is usually fast, too.

Comparison

The Efficient RANSAC algorithm is very quick, however, since it is not deterministic, some small shapes might be missed in the detection process.

Instead, the region growing algorithm usually takes longer to complete, but it may provide better quality output in the presence of large scenes with numerous small details. Since it iterates throughout all items of the scene, there are fewer chances to miss a shape. In addition, it is deterministic (for a given input and a given set of parameters, it always returns the same output, whereas the Efficient RANSAC algorithm is randomized and so the output may vary at each run), see Figure Figure 77.11.

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

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 exactly equal;
  • Parallel planes that are near coplanar are made exactly coplanar;
  • Planes that are near orthogonal are made exactly orthogonal;
  • Planes that are near symmetrical with respect to a user-defined axis are made exactly symmetrical.

The user can choose to regularize only one or several of these four 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 Shape_detection/efficient_RANSAC_and_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 <fstream>
#include <iostream>
#include <CGAL/property_map.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Point_with_normal_3.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Shape_detection/Efficient_RANSAC.h>
#include <CGAL/Regularization/regularize_planes.h>
// Type declarations.
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 ||
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.
Efficient_ransac efficient_ransac;
efficient_ransac.set_input(points);
efficient_ransac.add_shape_factory<Plane>();
efficient_ransac.detect();
Efficient_ransac::Plane_range planes = efficient_ransac.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;
}

Deprecated Components

The new version (see all examples above) of the class CGAL::Shape_detection::Region_growing is not compatible with the old API. For the old API, see an example below and use CGAL::Shape_detection::deprecated::Region_growing_depr instead.


File Shape_detection/shape_detection_basic_deprecated.cpp

#if defined (_MSC_VER) && !defined (_WIN64)
#pragma warning(disable:4244) // boost::number_distance::distance()
// converts 64 to 32 bits integers
#endif
#include <CGAL/internal/disable_deprecation_warnings_and_errors.h>
#include <fstream>
#include <iostream>
#include <CGAL/property_map.h>
#include <CGAL/IO/read_xyz_points.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Shape_detection_3.h>
// 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.
typedef CGAL::Shape_detection_3::Shape_detection_traits
<Kernel, Pwn_vector, Point_map, Normal_map> Traits;
typedef CGAL::Shape_detection_3::Efficient_RANSAC<Traits> Efficient_ransac;
typedef CGAL::Shape_detection_3::Region_growing_depr<Traits> Region_growing;
typedef CGAL::Shape_detection_3::Plane<Traits> Plane;
// This program both works for RANSAC and Region Growing.
// This example is using deprecated code!
// Please update your code to the new version using other examples!
template<typename ShapeDetection>
int run(const char* filename) {
// Points with normals.
Pwn_vector points;
// Load a 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::cout << "Error: cannot read the file cube.pwn" << std::endl;
return EXIT_FAILURE;
}
// Instantiate a shape detection engine.
ShapeDetection shape_detection;
// Provide the input data.
shape_detection.set_input(points);
// Register planar shapes via template method.
shape_detection.template add_shape_factory<Plane>();
// Detect registered shapes with default parameters.
shape_detection.detect();
// Print 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");
}

The old API is still supported, but will be removed in the next releases. In particular, the classes CGAL::Shape_detection::deprecated::Region_growing_depr and CGAL::Shape_detection::deprecated::Shape_detection_traits will be removed. Please update your code.

History

The Efficient RANSAC component was developed by Sven Oesau based on the prototype version created by Yannick Verdie, with the help of Clément Jamin and under the supervision of Pierre Alliez. Plane regularization was added by Simon Giraudot based on the prototype version developed by Florent Lafarge.

The region growing algorithm on a 3D point set was first implemented by Simon Giraudot based on the prototype version developed by Florent Lafarge and then generalized to arbitrary items including versions for a 2D point set, a 3D point set, and a polygon mesh by Thien Hoang during the Google Summer of Code 2018 under the supervision of Dmitry Anisimov.

Acknowledgments

The authors wish to thank our reviewers Andreas Fabri, Sébastien Loriot, and Simon Giraudot for helpful comments and discussions.