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

Amal Dev Parakkat, Stefan Ohrhallinger, Elmar Eisemann, Pooran Memari
The package provides a Delaunay triangulation based surface reconstruction algorithm from an unorganized point set. The reconstructed surface consists of the interface between Voronoi balls, which approximate the interior and exterior medial balls. This package contains the so-called global variant of the ball merge algorithm, which is carefully designed to reconstruct a watertight surface from real-world scanned data sets, exhibiting mild noise, outliers, and small gaps in the data.
Introduced in: CGAL 6.2
Depends on: 3D Triangulations
BibTeX: cgal:poem-bmsr-25b
License: GPL

Classified Reference Pages

Concepts

Functions

Class

Classes

class  CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag >
 Class that can be used for running the 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())
 creates two watertight meshes approximating the surface, and puts 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>

creates two watertight meshes approximating the surface, and puts the resulting triangle faces in out_triangles1 and out_triangles2.

As this function creates two shells (outer and inner) the input point set shall only be sampled on a single connected component. Output triangle faces are triple of indices refering to the position of the 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.