\( \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 4.8.2 - Poisson Surface Reconstruction
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
CGAL::Poisson_reconstruction_function< Gt > Class Template Reference

#include <CGAL/Poisson_reconstruction_function.h>

Definition

Implementation of the Poisson Surface Reconstruction method.

Given a set of 3D points with oriented normals sampled on the boundary of a 3D solid, the Poisson Surface Reconstruction method [2] solves for an approximate indicator function of the inferred solid, whose gradient best matches the input normals. The output scalar function, represented in an adaptive octree, is then iso-contoured using an adaptive marching cubes.

Poisson_reconstruction_function implements a variant of this algorithm which solves for a piecewise linear function on a 3D Delaunay triangulation instead of an adaptive octree.

Template Parameters
GtGeometric traits class.
Is Model Of:
ImplicitFunction
Examples:
Poisson_surface_reconstruction_3/poisson_reconstruction_example.cpp.

Types

typedef Gt Geom_traits
 Geometric traits class.
 
typedef Triangulation::Cell_handle Cell_handle
 
typedef Geom_traits::FT FT
 number type.
 
typedef Geom_traits::Point_3 Point
 point type.
 
typedef Geom_traits::Vector_3 Vector
 vector type.
 
typedef Geom_traits::Sphere_3 Sphere
 

Creation

template<typename InputIterator , typename PointPMap , typename NormalPMap >
 Poisson_reconstruction_function (InputIterator first, InputIterator beyond, PointPMap point_pmap, NormalPMap normal_pmap)
 Creates a Poisson implicit function from the range of points [first, beyond). More...
 

Operations

Sphere bounding_sphere () const
 Returns a sphere bounding the inferred surface.
 
template<class SparseLinearAlgebraTraits_d >
bool compute_implicit_function (SparseLinearAlgebraTraits_d solver, bool smoother_hole_filling=false)
 This function must be called after the insertion of oriented points. More...
 
FT operator() (const Point &p) const
 ImplicitFunction interface: evaluates the implicit function at a given 3D query point. More...
 
Point get_inner_point () const
 Returns a point located inside the inferred surface.
 

Constructor & Destructor Documentation

template<class Gt >
template<typename InputIterator , typename PointPMap , typename NormalPMap >
CGAL::Poisson_reconstruction_function< Gt >::Poisson_reconstruction_function ( InputIterator  first,
InputIterator  beyond,
PointPMap  point_pmap,
NormalPMap  normal_pmap 
)

Creates a Poisson implicit function from the range of points [first, beyond).

Template Parameters
InputIteratoriterator over input points.
PointPMapis a model of ReadablePropertyMap with a value_type = Point. It can be omitted if InputIterator value_type is convertible to Point.
NormalPMapis a model of ReadablePropertyMap with a value_type = Vector.
Parameters
firstiterator over the first input point.
beyondpast-the-end iterator over the input points.
point_pmapproperty map: value_type of InputIterator -> Point (the position of an input point).
normal_pmapproperty map: value_type of InputIterator -> Vector (the oriented normal of an input point).

Member Function Documentation

template<class Gt >
template<class SparseLinearAlgebraTraits_d >
bool CGAL::Poisson_reconstruction_function< Gt >::compute_implicit_function ( SparseLinearAlgebraTraits_d  solver,
bool  smoother_hole_filling = false 
)

This function must be called after the insertion of oriented points.

It computes the piecewise linear scalar function operator() by: applying Delaunay refinement, solving for operator() at each vertex of the triangulation with a sparse linear solver, and shifting and orienting operator() such that it is 0 at all input points and negative inside the inferred surface.

Template Parameters
SparseLinearAlgebraTraits_dSymmetric definite positive sparse linear solver. If Eigen 3.1 (or greater) is available and CGAL_EIGEN3_ENABLED is defined, an overload with Eigen_solver_traits<Eigen::ConjugateGradient<Eigen_sparse_symmetric_matrix<double>::EigenType> > as default solver is provided.
Parameters
solversparse linear solver.
smoother_hole_fillingcontrols if the Delaunay refinement is done for the input points, or for an approximation of the surface obtained from a first pass of the algorithm on a sample of the points.
Returns
false if the linear solver fails.
template<class Gt >
FT CGAL::Poisson_reconstruction_function< Gt >::operator() ( const Point p) const

ImplicitFunction interface: evaluates the implicit function at a given 3D query point.

The function compute_implicit_function() must be called before the first call to operator().