- the value of \( \delta \)
- Type: double
- Default: 1.7
|
CGAL 6.0 - Ball Merge Surface Reconstruction
|
#include <CGAL/Ball_merge_surface_reconstruction.h>
Class that can be used for running the ball merge surface reconstruction algorithms.
Once input points passed, it is possible to run a reconstruction algorithm with different parameters without rebuilding the internal Delaunay Triangulation.
| Traits | a model of DelaunayTriangulationTraits_3, with default using the value type of PointRange plugged in Kernel_traits |
| ConcurrencyTag | enables sequential versus parallel algorithm. Possible values are Sequential_tag, Parallel_tag, and Parallel_if_available_tag. |
Public Member Functions | |
| template<class PointRange > | |
| void | build_triangulation (const PointRange &points) |
| sets the input points for the triangulation, and (build the internal triangulation. | |
| template<class NamedParameters = parameters::Default_named_parameters, class TripleIndexRange > | |
| void | local_reconstruction (TripleIndexRange &out_triangles, const NamedParameters &np=parameters::default_values()) |
creates a triangle soup approximating the surface with sample points passed to build_triangulation(), and puts the resulting triangule faces in out_triangles. | |
| template<class NamedParameters = parameters::Default_named_parameters, class TripleIndexRange > | |
| void | global_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 triangule faces in out_triangles1 and out_triangles2. | |
| void CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag >::build_triangulation | ( | const PointRange & | points | ) |
sets the input points for the triangulation, and (build the internal triangulation.
If called several times, only the last point range will be considered for the reconstructions.
| PointRange | a model of RandomAccessContainer, with Traits::Point_3 as value type |
| points | is the input points |
| void CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag >::global_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 triangule faces in out_triangles1 and out_triangles2.
| TripleIndexRange | a model of BackInsertionSequence with value_type being a model of RandomAccessContainer and BackInsertionSequence with value_type being constructible from std::size_t. |
| out_triangles1 | is the output parameter storing the first resulting mesh |
| out_triangles2 | is the output parameter storing the second resulting mesh |
| np | an optional sequence of Named Parameters among the ones listed below |
|
| void CGAL::Ball_merge_surface_reconstruction< Traits, ConcurrencyTag >::local_reconstruction | ( | TripleIndexRange & | out_triangles, |
| const NamedParameters & | np = parameters::default_values() |
||
| ) |
creates a triangle soup approximating the surface with sample points passed to build_triangulation(), and puts the resulting triangule faces in out_triangles.
| TripleIndexRange | a model of BackInsertionSequence with value_type being a model of RandomAccessContainer and BackInsertionSequence with value_type being constructible from std::size_t. |
| out_triangles | is the output parameter storing triangles approximating the surface |
| np | an optional sequence of Named Parameters among the ones listed below |
| |
|