CGAL 6.1 - Ball Merge Surface Reconstruction
Loading...
Searching...
No Matches
CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag > Class Template Reference

#include <CGAL/Ball_merge_surface_reconstruction.h>

Definition

template<typename Traits, typename ConcurrencyTag>
class CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag >

Class that can be used for running the ball merge surface reconstruction algorithms.

Once input points are passed, it is possible to run a reconstruction algorithm with different parameters without rebuilding the internal Delaunay triangulation.

Template Parameters
Traitsa model of DelaunayTriangulationTraits_3, with default using the value type of PointRange plugged in Kernel_traits
ConcurrencyTagenables sequential versus parallel algorithm. Possible values are Sequential_tag, Parallel_tag, and Parallel_if_available_tag.

Public Member Functions

 Ball_merge_surface_reconstruction ()
 default constructor
 
template<class PointRange >
 Ball_merge_surface_reconstruction (const PointRange &points)
 constructs the class and calls build_triangulation(points) .
 
template<class PointRange >
void build_triangulation (const PointRange &points)
 sets the input points for the triangulation, and builds the internal triangulation.
 
template<class NamedParameters = parameters::Default_named_parameters, class TripleIndexRange >
double reconstruction (TripleIndexRange &out_triangles1, TripleIndexRange &out_triangles2, const NamedParameters &np=parameters::default_values())
 creates two watertight meshes approximating the surface with sample points passed to build_triangulation(), and puts the resulting triangle faces in out_triangles1 and out_triangles2.
 

Constructor & Destructor Documentation

◆ Ball_merge_surface_reconstruction()

template<typename Traits , typename ConcurrencyTag >
template<class PointRange >
CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag >::Ball_merge_surface_reconstruction ( const PointRange &  points)

constructs the class and calls build_triangulation(points) .

Template Parameters
PointRangea model of RandomAccessContainer, with Traits::Point_3 as value type
Parameters
pointsis the input points

Member Function Documentation

◆ build_triangulation()

template<typename Traits , typename ConcurrencyTag >
template<class PointRange >
void CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag >::build_triangulation ( const PointRange &  points)

sets the input points for the triangulation, and builds the internal triangulation.

If called several times, only the last point range will be considered for the reconstructions.

Template Parameters
PointRangea model of RandomAccessContainer, with Traits::Point_3 as value type
Parameters
pointsis the input points

◆ reconstruction()

template<typename Traits , typename ConcurrencyTag >
template<class NamedParameters = parameters::Default_named_parameters, class TripleIndexRange >
double CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag >::reconstruction ( TripleIndexRange &  out_triangles1,
TripleIndexRange &  out_triangles2,
const NamedParameters &  np = parameters::default_values() 
)

creates two watertight meshes approximating the surface with sample points passed to build_triangulation(), and puts the resulting triangle faces in out_triangles1 and out_triangles2.

Output triangle faces are triple of indices refering to the position of the input points passed to build_triangulation().

Note
As this function creates two shells (outer and inner in arbitrary order) the input point set must only be sampled on a single connected component.
Template Parameters
TripleIndexRangea model of BackInsertionSequence with value_type being a model of RandomAccessContainer and BackInsertionSequence with value_type being constructible from std::size_t.
Parameters
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

Precondition
build_triangulation() should have been called before calling this function.