CGAL 6.2 - Ball Merge Surface Reconstruction
Loading...
Searching...
No Matches
Ball Merge Surface Reconstruction Reference

Amal Dev Parakkat, Stefan Ohrhallinger, Elmar Eisemann, Pooran Memari
This package offers a surface reconstruction algorithm based on Delaunay triangulation, designed to work with unstructured point clouds. It builds a watertight surface by merging balls centered at Voronoi vertices, which approximate the medial axisof the shape. The algorithm uses a global variant of the ball merge strategy, making it robust to typical issues in scanned data—such as noise, outliers, and small gaps. It’s especially suited for reconstructing clean, closed surfaces from real-world 3D scans.
Introduced in: CGAL 6.2
Depends on: 3D Triangulations
BibTeX: cgal:poem-bmsr-26a
License: GPL

Classified Reference Pages

Concepts

Functions

Class

Classes

class  CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag >
 This class provides an interface for executing ball merge surface reconstruction algorithms. More...
 

Functions

template<class NamedParameters = parameters::Default_named_parameters, class PointRange , class TripleIndexRange >
double CGAL::ball_merge_surface_reconstruction (const PointRange &points, TripleIndexRange &out_triangles1, TripleIndexRange &out_triangles2, const NamedParameters &np=parameters::default_values())
 generates two meshes approximating the input surface and stores the resulting triangle faces in out_triangles1 and out_triangles2
 

Function Documentation

◆ ball_merge_surface_reconstruction()

template<class NamedParameters = parameters::Default_named_parameters, class PointRange , class TripleIndexRange >
double CGAL::ball_merge_surface_reconstruction ( const PointRange &  points,
TripleIndexRange &  out_triangles1,
TripleIndexRange &  out_triangles2,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Ball_merge_surface_reconstruction.h>

generates two meshes approximating the input surface and stores the resulting triangle faces in out_triangles1 and out_triangles2

Precondition

This function constructs two shells (outer and inner). Therefore, the input point set must sample a single connected component.

Output Format

Each triangle face is represented as a triplet of indices referring to the positions of the corresponding input points in points.

Template Parameters
PointRangea model of the concepts RandomAccessContainer
TripleIndexRangea model of BackInsertionSequence with value_type being a model of RandomAccessContainer and BackInsertionSequence with value_type being constructible from std::size_t.
Parameters
pointsis the input points representing a single connected component of a watertight surface
out_triangles1is the output parameter storing the first resulting mesh
out_triangles2is the output parameter storing the second resulting mesh
npan optional sequence of Named Parameters among the ones listed below
Returns
the delta value used for the reconstruction
Optional Named Parameters
  • the value of \( \delta \)
  • Type: double
  • Default: Value iteratively estimated, as described in Parameter Settings
Examples
Ball_merge_surface_reconstruction/ball_merge_reconstruction.cpp.